REQUIREMENT: Before proceeding, install or update your stenv distribution. stenv is the replacement for AstroConda, which is unsupported as of February 2023.

This notebook currently fails to execute, use as reference only

ACS Image Reduction for Subarray Data#

Introduction#


Subarray data requires different considerations than working with full frame images. This notebook will guide you through working with standard post-SM4 subarray data and custom subarray data.

After Servicing Mission 4 (SM4; May 2009), the installation of the ASIC during the repair of ACS introduced \(1/f\) noise in all ACS images. In the calacs pipeline, only full-frame data have this striping removed. To correct subarray data, the alternative acs_destripe_plus pipeline must be used, which will apply all of calibration steps normally performed by calacs in addition to de-striping for subarray data. De-striping is only possible for 2K subarrays after SM4 until Cycle 24, after which a change to the flight software makes all subarrays eligible for de-striping.

This tutorial will show you how to handle…#

1. Post-SM4 Subarray Data#

  • Update header keywords.

  • Clean Subarray images with the option to correct CTE losses.

  • Update WCS solution.

2. Custom Subarray Data#

  • Use AstroDrizzle with ASN files to combine images.

Imports#


Here we list the Python packages used in this notebook. Links to the documentation for each module is provided for convenience.

Package Name

module

docs

used for

os

system

link

command line input

os

environ

link

setting environments

shutil

rmtree

link

remove directory tree

glob

glob

link

search for files based on Unix shell rules

astroquery.mast

Observations

link

download data from MAST

astropy.io

fits

link

access and update fits files

astropy.table

Table

link

constructing and editing in a tabular format

stwcs

updatewcs

link

update wcs solution

acstools

acs_destripe_plus

link

destripe acs images and optionally CTE correct

acstools

utils_calib

link

check for settings in oscntab

import os
import shutil

from astroquery.mast import Observations
from astropy.io import fits
from astropy.table import Table

from stwcs import updatewcs
from acstools import (acs_destripe_plus, utils_calib)

Here we set environment variables for later use with the Calibration Reference Data System (CRDS).

os.environ['CRDS_SERVER_URL'] = 'https://hst-crds.stsci.edu'
os.environ['CRDS_SERVER'] = 'https://hst-crds.stsci.edu'
os.environ['CRDS_PATH'] = './crds_cache'
os.environ['jref'] = './crds_cache/references/hst/acs/'

Download the Data#


Here we download all of the data required for this notebook. This is an important step! Some of the image processing steps require all relevant files to be in the working directory. We recommend working with a brand new directory for every new set of data.

GO Proposal 14511: “ACS CCD Stability Monitor”#

Our data for the first example comes from a routine calibration program of two images using the post-Cycle 24 WFC1A-1K 1024 x 2072 subarray. We will only use the data associated with the observation id JD5702JWQ.

Using the python package astroquery, we can download files from the MAST archive.

MAY CHANGE: The argument "mrp_only" stands for "minimum recommended products only". It currently needs to be set to False, although in the future, False is intended to be set as the default and can be left out.
obs_table = Observations.query_criteria(proposal_id=14511, obs_id='JD5702JWQ')

dl_table = Observations.download_products(obs_table['obsid'],
                                          productSubGroupDescription=['RAW'],
                                          mrp_only=False)

We’ll use the package os to put all of these files in our working directory for convenience.

for row in dl_table:
    oldfname = row['Local Path']
    newfname = os.path.basename(oldfname)
    os.rename(oldfname, newfname)

GO Proposal 10206: “What drives the outflows in powerful radio galaxies?”#

For the second example, we will use a ramp filter observation of the galaxy PLS1345+12 (HST proposal 10206). The association name is J92SA0010, and we will only use one image in the association: J92SA0W6Q.

Again, we use astroquery, to download files from the MAST archive. Along with the raw images, we will also need the spt (telemetry) files to reduce our custom subarray data.

obs_table = Observations.query_criteria(proposal_id=10206, obs_id='J92SA0010')

dl_table = Observations.download_products(obs_table['obsid'],
                                          productSubGroupDescription=['RAW', 'SPT'],
                                          mrp_only=False)

Again, we’ll use the package os to put all of these files in our working directory for convenience.

for row in dl_table:
    oldfname = row['Local Path']
    newfname = os.path.basename(oldfname)
    os.rename(oldfname, newfname)

Now that all of our files are in the current working directory, we delete the leftover MAST file structure.

shutil.rmtree('mastDownload')

File Information#

The structure of the fits files from ACS may be different depending on what kind of observation was made. For more information, refer to Section 2.2 of the ACS Data Handbook.

Raw Files#

A standard raw image file from a subarray has the same structure as you’d expect from full frame observation from ACS/WCS.

Ext

Name

Type

Contains

0

PRIMARY

(PrimaryHDU)

Meta-data related to the entire file.

1

SCI (Image)

(ImageHDU)

Raw image data.

2

ERR (Error)

(ImageHDU)

Error array.

3

DQ (Data Quality)

(ImageHDU)

Data quality array.

SPT Files#

SPT files contain telemetry and engineering data from the telescope.

Ext

Name

Type

Contains

0

PRIMARY

(PrimaryHDU)

Meta-data related to the entire file.

1

UDL (Image)

(ImageHDU)

Raw image data.

You can always use .info() on an HDUlist for an overview of the structure

with fits.open('j92sa0w6q_raw.fits') as hdulist:
    hdulist.info()
    
with fits.open('j92sa0w6q_spt.fits') as hdulist:
    hdulist.info()

De-striping and CTE Correction of Post-SM4 Subarray Observations#


Download Calibration Files: Update Header Keywords#

We can call the Calibration Reference Data System (CRDS) to get the associated calibration files for this image.

First we will need to turn on the CTE correction switch in the primary header. Turning on this switch will notify the CRDS bestrefs tool to add the names of the CTE correction parameters table ‘PCTETAB’ and CTE-corrected dark current image ‘DRKCFILE’ reference files to the header.

sbc_fits = 'jd5702jwq_raw.fits'

with fits.open(sbc_fits, mode='update') as hdulist:
    hdulist[0].header['PCTECORR'] = 'PERFORM'

Download Calibration Files#

The CRDS program can be run from the terminal to download the files specific to your observation. This program is packaged with astroconda.

The command line input to download our files for jd5702jwq_raw.fits is as follows:

crds bestrefs --files jd5702jwq_raw.fits --sync-references=1 --update-bestrefs

Here, we use the python package os to run this command.

cmd_input = 'crds bestrefs --files {:} --sync-references=1 --update-bestrefs'.format(sbc_fits)
os.system(cmd_input)

Clean Subarray Images#

Next we will run the acs_destripe_plus code on our image. This will execute all of the calibration steps that are set to ‘PERFORM’ in the primary header of the FITS file.

The acs_destripe_plus code will produce the FLT file is the calibrated output product from CALACS, jd5702jmq_flt.fits. With the CTE correction turned on, an FLC file will also be produced, which is the same as the FLT file but with the CTE correction applied.

acs_destripe_plus.destripe_plus(sbc_fits, cte_correct=True)

Correct the WCS#

The subarray products produced by this process do not have the proper WCS information in the header. The WCS is normally updated by the pipeline via an additional call to AstroDrizzle. Here, we can manually update the WCS of our FLC product using stwcs.updatewcs.

updatewcs.updatewcs('jd5702jwq_flc.fits', use_db=False)

Reducing Custom Subarray Data#


Download Calibration Files#

Like before, we can use CRDS to get the associated calibration files for this image.

sbc_fits = 'j92sa0w6q_raw.fits'
cmd_input = 'crds bestrefs --files {:} --sync-references=1 --update-bestrefs'.format(sbc_fits)
os.system(cmd_input)

Access OSCN Table#

The name in the header is of the format ‘ref\\(oscn_name.fits</font>', therefore we need to split the string on the '\)’ character.

prihdr = fits.getheader(sbc_fits)
scihdr = fits.getheader(sbc_fits, ext=1)

We also want to retrieve the OSCNTAB reference file from the JREF directory. We can get the name of the file from the primary header of the image.

oscn_name = prihdr['OSCNTAB'].split('$')[-1]
path_to_oscn = os.path.join(os.environ['jref'], oscn_name)

print(path_to_oscn)

The utils_calib.check_oscntab from acstools checks the OSCNTAB file if any entry matches the combination of readout amplifier, image size, and number of bias prescan columns for a given subarray observation. We will need to use several header keyword values to test if a subarray is in OSCNTAB.

oscnrec = fits.getdata(path_to_oscn)
oscnhdr = fits.getheader(path_to_oscn)

oscntable = Table(oscnrec)

The raw image has 1180 columns and 1200 rows, which does not correspond to any entry in the OSCNTAB file, but a visual examination of the image shows that it contains bias prescan columns.

Get the Bias Prescan Information#

From the science image PRI and SCI extension headers that we opened earlier, we can get the information about the readout amplifier and dimensions of the image.

amp = prihdr['CCDAMP']
xsize = scihdr['NAXIS1']
ysize = scihdr['NAXIS2']

print(xsize, ysize)

To get information on the number of prescan columns (if any), we need to access the SPT first extension header.

spthdr = fits.getheader('j92sa0w6q_spt.fits', ext=1)

leading = spthdr['OVERSCNA']
trailing = spthdr['OVERSCNB']

Finally, we check if this subarray definition is in the OSCNTAB file. The code returns a boolean result, which we have saved as the variable supported, to describe this.

supported = utils_calib.check_oscntab(path_to_oscn, amp, xsize, ysize, leading, trailing)
print(supported)

Update OSCNTAB#

Now that we have confirmed that the OSCNTAB file does not contain information about our subarray data, we need to add a new row to the table with our definitions. Let’s first view the first few rows of OSCNTAB to see what our new entry needs to look like.

print(oscntable)

We can also choose to just view all of the columns names

print(oscntable.colnames)

Several column names are obvious, but here we define the less obvious ones.

Column Name

Description

BINX, BINY

Binning in X and Y. ACS data are never binned, so these will always be 1.

TRIMXn

Number of prescan columns on the left (1) and right (2) sides of the image to remove.

TRIMYn

Number of virtual rows on the bottom (1) and top (2) sides of the image to remove. For subarray data, these are always 0.

BIASSECTAn

Start and end columns to use for the bias level estimation on the left (A) side of the image.

BIASSECTBn

Start and end columns to use for the bias level estimation on the right (B) side of the image.

VXn, VYn

The coordinates of the bottom-left (VX1, VY1) and top-right (VX2, VY2) corners of the virtual overscan region

The following line sets chip to 1 if the subarray is on WFC1, and 2 if the subarray is on WFC2.

chip = 1 if amp in ['A', 'B'] else 2

For the BIASSECTAn and BIASSECTBn values, we want to use the six columns of prescan nearest the exposed area of the CCD.

bias_a = [0, 0]
if leading > 0:
    bias_a[1] = leading
    bias_a[0] = leading-5 if leading > 5 else 0
        
bias_b = [0, 0]
if trailing > 0:
    bias_b[0] = xsize+1
    bias_b[1] = xsize+5 if trailing > 5 else xsize+trailing

Now we can define a new row for our settings. For subarray data, as there is no virtual overscan, the VXn and VYn values will always be 0. Here we use a dictionary to explicitly define the new values.

new_row = {'CCDAMP': amp,
           'CCDCHIP': chip, 
           'BINX': 1, 
           'BINY': 1, 
           'NX': xsize,
           'NY': ysize, 
           'TRIMX1': leading, 
           'TRIMX2': trailing, 
           'TRIMY1': 0, 
           'TRIMY2': 0, 
           'BIASSECTA1': bias_a[0], 
           'BIASSECTA2': bias_a[1], 
           'BIASSECTB1': bias_b[0], 
           'BIASSECTB2': bias_b[1], 
           'DESCRIPTION': 'Custom OSCN', 
           'VX1': 0, 
           'VX2': 0, 
           'VY1': 0, 
           'VY2': 0}

Now we need to open the custom OSCTNAB file and update the table to have our new definition. We will also need to update the raw FITS file to point to the new custom OSCNTAB reference file.

oscntable.add_row(new_row)
print(oscntable)

Now we have an idential FITS rec table, but with an additional row for our new information.

oscnrec_new = fits.FITS_rec.from_columns(oscnrec.columns, 
                                         nrows=len(oscntable))
print(Table(oscnrec_new))

Let’s populate that last row with our new data!

oscnrec_new[-1] = tuple(oscntable[-1])

Reprint the table and check that we have entered it correctly.

print(Table(oscnrec_new))

Now we can open up our oscn table and replace the old fits rec with our new one

with fits.open(path_to_oscn, mode='update') as hdu:
    hdu[1].data = oscnrec_new

Using the same check we did earlier, we can use check_oscntab to see if our settings are defined. If everything was done correctly, the following line should print “True”!

supported = utils_calib.check_oscntab(path_to_oscn, amp, xsize, ysize, leading, trailing)

print('Defined in OSCNTAB? {}'.format(supported))

For more help:#

More details may be found on the ACS website and in the ACS Instrument and Data Handbooks.

Please visit the HST Help Desk. Through the help desk portal, you can explore the HST Knowledge Base and request additional help from experts.


Top of Page Space Telescope Logo