GrismFLT¶
- class grizli.model.GrismFLT(grism_file='', sci_extn=1, direct_file='', pad=(64, 256), ref_file=None, ref_ext=0, seg_file=None, shrink_segimage=True, force_grism='G141', verbose=True, process_jwst_header=True, use_jwst_crds=False)[source]¶
Bases:
objectScripts for modeling of individual grism FLT images
Read FLT files and, optionally, reference/segmentation images.
- Parameters
- grism_filestr
Grism image (optional). Empty string or filename of a FITS file that must contain extensions (‘SCI’,
sci_extn), (‘ERR’,sci_extn), and (‘DQ’,sci_extn). For example, a WFC3/IR “FLT” FITS file.- sci_extnint
EXTNAME of the file to consider. For WFC3/IR this can only be 1. For ACS and WFC3/UVIS, this can be 1 or 2 to specify the two chips.
- direct_filestr
Direct image (optional). Empty string or filename of a FITS file that must contain extensions (‘SCI’,
sci_extn), (‘ERR’,sci_extn), and (‘DQ’,sci_extn). For example, a WFC3/IR “FLT” FITS file.- padint, int
Padding to add around the periphery of the images to allow modeling of dispersed spectra for objects that could otherwise fall off of the direct image itself. Modeling them requires an external reference image (
ref_file) that covers an area larger than the individual direct image itself (e.g., a mosaic of a survey field).For WFC3/IR spectra, the first order spectra reach 248 and 195 pixels for G102 and G141, respectively, and
padcould be set accordingly if the reference image is large enough.- ref_filestr or
ImageHDU/PrimaryHDU Image mosaic to use as the reference image in place of the direct image itself. For example, this could be the deeper image drizzled from all direct images taken within a single visit or it could be a much deeper/wider image taken separately in perhaps even a different filter.
Note
Assumes that the WCS are aligned between
grism_file,direct_fileandref_file!- ref_extint
FITS extension to use if
ref_fileis a filename string.- seg_filestr or
ImageHDU/PrimaryHDU Segmentation image mosaic to associate pixels with discrete objects. This would typically be generated from a rectified image like
ref_file, though here it is not required thatref_fileandseg_filehave the same image dimensions but rather just that the WCS are aligned between them.- shrink_segimagebool
Try to make a smaller cutout of the reference images to speed up blotting and array copying. This is most helpful for very large input mosaics.
- force_grismstr
Use this grism in “simulation mode” where only
direct_fileis specified.- verbosebool
Print status messages to the terminal.
- process_jwst_headerbool
Process JWST headers to get the grism configuration file.
- use_jwst_crdsbool
Use the
jwst.datamodelspackage to read JWST headers and configuration files. If False, use thegrizli.grismconfpackage.
- Attributes
- grism, direct
ImageData Grism and direct image data and parameters
- conf
aXeConf Grism configuration object.
- segarray-like
Segmentation image array.
- modelarray-like
Model of the grism exposure with the same dimensions as the full detector array.
- object_dispersersdict
Container for storing information about what objects have been added to the model of the grism exposure
- catalog
Table Associated photometric catalog. Not required.
- grism, direct
Methods Summary
apply_POM([warn_if_too_small, verbose])Apply pickoff mirror mask to segmentation map to control sources that are dispersed onto the detector
blot_catalog(input_catalog[, columns, ...])Compute detector-frame coordinates of sky positions in a catalog.
compute_full_model([ids, mags, mag_limit, ...])Compute flat-spectrum model for multiple objects.
compute_model_orders([id, x, y, size, mag, ...])Compute dispersed spectrum for a given object id
get_dispersion_PA([decimals])Compute exact PA of the dispersion axis, including tilt of the trace and the FLT WCS
get_trace_region_from_sky(ra, dec[, width])Make a region file for the trace in pixel coordinates given sky position TBD
load_from_fits(save_file)Load saved data from a FITS file
load_photutils_detection([seg_file, ...])Load segmentation image and catalog, either from photutils or SExtractor.
mask_mosaic_edges([sky_poly, verbose, ...])Mask edges of exposures that might not have modeled spectra
old_make_edge_mask([scale, force])Make a mask for the edge of the grism FoV that isn't covered by the direct image
photutils_detection([use_seg, data_ext, ...])Use photutils to detect objects and make segmentation map
process_ref_file(ref_file[, ref_ext, ...])Read and blot a reference image
process_seg_file(seg_file[, ...])Read and blot a rectified segmentation image
save_full_pickle([verbose])Save entire
GrismFLTobject to a picklesave_model([overwrite, verbose])Save model properties to FITS file
save_wcs([overwrite, verbose])Save WCS of the direct and grism images to separate FITS files.
smooth_mask([gaussian_width, threshold])Compute a mask where smoothed residuals greater than some value
transform_JWST_WFSS([verbose])Rotate data & wcs so that spectra are increasing to +x
Methods Documentation
- apply_POM(warn_if_too_small=True, verbose=True)[source]¶
Apply pickoff mirror mask to segmentation map to control sources that are dispersed onto the detector
- Parameters
- warn_if_too_smallbool
Print a warning if the
padparameter is too small to cover the full POM region.- verbosebool
Print status messages to the terminal.
- Returns
- True if completed successfully.
- blot_catalog(input_catalog, columns=['id', 'ra', 'dec'], sextractor=False, ds9=None)[source]¶
Compute detector-frame coordinates of sky positions in a catalog.
- Parameters
- input_catalog
Table Full catalog with sky coordinates. Can be SExtractor or other.
- columns[str,str,str]
List of columns that specify the object id, R.A. and Decl. For catalogs created with SExtractor this might be [‘NUMBER’, ‘X_WORLD’, ‘Y_WORLD’].
Detector coordinates will be computed with
self.direct.wcs.all_world2pixwithorigin=1.- sextractorbool
If True, then the input catalog is assumed to be from sextractor and the columns are [‘NUMBER’, ‘X_WORLD’, ‘Y_WORLD’].
- ds9
DS9, optional If provided, load circular regions at the derived detector coordinates.
- input_catalog
- Returns
- catalog
Table New catalog with columns ‘x_flt’ and ‘y_flt’ of the detector coordinates. Also will copy the
columnsnames to columns with names ‘id’,’ra’, and ‘dec’ if necessary, e.g., for SExtractor catalogs.
- catalog
- compute_full_model(ids=None, mags=None, mag_limit=22, store=True, verbose=False, size=10, min_size=26, compute_size=True)[source]¶
Compute flat-spectrum model for multiple objects.
- Parameters
- idsNone, list, or
array id numbers to compute in the model. If None then take all ids from unique values in
self.seg.- magsNone, float, or list /
array magnitudes corresponding to list if
ids. If None, then compute magnitudes based on the flux in segmentation regions and zeropoints determined from PHOTFLAM and PHOTPLAM.- mag_limitfloat
Limiting magnitude for computing the model. If
magsis None, then compute magnitudes for all objects in the segmentation image and only compute the model for objects withmagsless thanmag_limit.- storebool
Store the computed models in the
self.object_dispersersattribute.- verbosebool
Print status messages to the terminal.
- size, compute_sizeint, bool
Sizes of individual cutouts, see
compute_model_orders.- min_sizeint
Minimum size of the cutout to extract. This is enforced when
compute_sizeis True.
- idsNone, list, or
- Returns
- Updated model stored in
self.modelattribute.
- Updated model stored in
- compute_model_orders(id=0, x=None, y=None, size=10, mag=-1, spectrum_1d=None, is_cgs=False, compute_size=False, max_size=None, min_size=26, store=True, in_place=True, get_beams=None, psf_params=None, verbose=True)[source]¶
Compute dispersed spectrum for a given object id
- Parameters
- idint
Object ID number to match in the segmentation image
- x, yfloat
Center of the cutout to extract
- sizeint
Radius of the cutout to extract. The cutout is equivalent to
>>> xc, yc = int(x), int(y) >>> thumb = self.direct.data['SCI'][yc-size:yc+size, xc-size:xc+size]
- magfloat
Specified object magnitude, which will be compared to the “MMAG_EXTRACT_[BEAM]” parameters in
self.confto decide if the object is bright enough to compute the higher spectral orders. Default of -1 means compute all orders listed inself.conf.beams- spectrum_1dNone or [
array,array] Template 1D spectrum to convolve with the grism disperser. If None, assumes trivial spectrum flat in f_lambda flux densities. Otherwise, the template is taken to be
>>> wavelength, flux = spectrum_1d
- is_cgsbool
Flux units of
spectrum_1d[1]are cgs f_lambda flux densities, rather than normalized in the detection band.- compute_sizebool
Ignore
x,y, andsizeand compute the extent of the segmentation polygon directly usingcompute_segmentation_limits.- max_sizeint or None
Enforce a maximum size of the cutout when using
compute_size.- min_sizeint
Minimum size of the cutout to extract. This is enforced when
compute_sizeis True.- storebool
If True, then store the computed beams in the OrderedDict
self.object_dispersers[id].If many objects are computed, this can be memory intensive. To save memory, set to False and then the function just stores the input template spectrum (
spectrum_1d) and the beams will have to be recomputed if necessary.- in_placebool
If True, add the computed spectral orders into
self.model. Otherwise, make a clean array with only the orders of the given object.- get_beamslist or None
Spectral orders to retrieve with names as defined in the configuration files, e.g., [‘A’] generally for the +1st order of HST grisms. If
None, then get all orders listed in thebeamsattribute of theaXeConfconfiguration object.- psf_paramslist
Optional parameters for generating an
EffectivePSFobject for the spatial morphology.- verbosebool
Print status messages to the terminal.
- Returns
- outputbool or
numpy.array If
in_placeis True, return status of True if everything goes OK. The computed spectral orders are stored in place inself.model.Returns False if the specified
idis not found in the segmentation array independent ofin_place.If
in_placeis False, return a full array including the model for the single object.
- outputbool or
- get_dispersion_PA(decimals=0)[source]¶
Compute exact PA of the dispersion axis, including tilt of the trace and the FLT WCS
- Parameters
- decimalsint or None
Number of decimal places to round to, passed to
round. If None, then don’t round.
- Returns
- dispersion_PAfloat
PA (angle East of North) of the dispersion axis.
- get_trace_region_from_sky(ra, dec, width=2)[source]¶
Make a region file for the trace in pixel coordinates given sky position TBD
- load_from_fits(save_file)[source]¶
Load saved data from a FITS file
- Parameters
- save_filestr
Filename of the saved output
- Returns
- True if completed successfully
- load_photutils_detection(seg_file=None, seg_cat=None, catalog_format='ascii.commented_header')[source]¶
Load segmentation image and catalog, either from photutils or SExtractor.
If SExtractor, use
catalog_format='ascii.sextractor'.- Parameters
- seg_filestr
Filename of the segmentation image.
- seg_catstr
Filename of the segmentation catalog.
- catalog_formatstr
Format of the catalog file.
- mask_mosaic_edges(sky_poly=None, verbose=True, force=False, err_scale=10, dq_mask=False, dq_value=1024, resid_sn=7)[source]¶
Mask edges of exposures that might not have modeled spectra
- Parameters
- sky_poly
Polygon Polygon in sky coordinates that defines the region to mask.
- verbosebool
Print status messages to the terminal.
- forcebool
Force the mask to be applied even if it has already been set.
- err_scalefloat
Scale factor to multiply to the mask before it’s applied to the
self.grism.data['ERR']array.- dq_maskbool
Mask the DQ array instead of the ERR array.
- dq_valueint
DQ value to apply to the masked pixels.
- resid_snfloat
Mask pixels where the residual is greater than this value times the error array.
- sky_poly
- old_make_edge_mask(scale=3, force=False)[source]¶
Make a mask for the edge of the grism FoV that isn’t covered by the direct image
- Parameters
- scalefloat
Scale factor to multiply to the mask before it’s applied to the
self.grism.data['ERR']array.- forcebool
Force apply the mask even if
self.has_edge_maskis set indicating that the function has already been run.
- Returns
- Nothing, updates
self.grism.data['ERR']in place. - Sets
self.has_edge_mask = True.
- Nothing, updates
- photutils_detection(use_seg=False, data_ext='SCI', detect_thresh=2.0, grow_seg=5, gauss_fwhm=2.0, verbose=True, save_detection=False, ZP=None)[source]¶
Use photutils to detect objects and make segmentation map
- Parameters
- use_segbool
Use the segmentation image in
self.segas the starting point for the detection image.- data_extstr
Extension of the
self.direct.dataarray to use for detection.- detect_threshfloat
Detection threshold, in sigma
- grow_segint
Number of pixels to grow around the perimeter of detected objects witha maximum filter
- gauss_fwhmfloat
FWHM of Gaussian convolution kernel that smoothes the detection image.
- verbosebool
Print logging information to the terminal
- save_detectionbool
Save the detection images and catalogs
- ZPfloat or None
AB magnitude zeropoint of the science array. If
Nonethen, try to compute based on PHOTFLAM and PHOTPLAM values and use zero if that fails.
- Returns
- statusbool
True if completed successfully. False if
data_ext=='REF'but no reference image found.- Stores an astropy.table.Table object to
self.catalogand a - segmentation array to
self.seg.
- process_ref_file(ref_file, ref_ext=0, shrink_segimage=True, verbose=True)[source]¶
Read and blot a reference image
- Parameters
- ref_filestr or
ImageHDU/PrimaryHDU Filename or
astropy.io.fitsImage HDU of the reference image.- ref_extint
FITS extension to use if
ref_fileis a filename string.- shrink_segimagebool
Try to make a smaller cutout of the reference image to speed up blotting and array copying. This is most helpful for very large input mosaics.
- verbosebool
Print some status information to the terminal
- ref_filestr or
- Returns
- statusbool
False if
ref_fileis None. True if completes successfully.- The blotted reference image is stored in the array attribute
self.direct.data['REF'].- The
ref_filterattribute is determined from the image header and the ref_photflamscaling is taken either from the header if possible, or- the global
photflamvariable defined at the top of this file.
- process_seg_file(seg_file, shrink_segimage=True, verbose=True)[source]¶
Read and blot a rectified segmentation image
- Parameters
- seg_filestr or
ImageHDU/PrimaryHDU Filename or
astropy.io.fitsImage HDU of the segmentation image.- shrink_segimagebool
Try to make a smaller cutout of the segmentation image to speed up blotting and array copying. This is most helpful for very large input mosaics.
- verbosebool
Print some status information to the terminal
- seg_filestr or
- Returns
- The blotted segmentation image is stored in the attribute
GrismFLT.seg.
- The blotted segmentation image is stored in the attribute
- save_model(overwrite=True, verbose=True)[source]¶
Save model properties to FITS file
- Parameters
- overwritebool
Overwrite existing files.
- verbosebool
Print status messages to the terminal.
- save_wcs(overwrite=True, verbose=True)[source]¶
Save WCS of the direct and grism images to separate FITS files.
- Parameters
- overwritebool
Overwrite existing files.
- verbosebool
Print status messages to the terminal.
- Returns
- Nothing, but writes WCS to files like
j123456m123456.01.FLT.wcs.fits.
- Nothing, but writes WCS to files like
- smooth_mask(gaussian_width=4, threshold=2.5)[source]¶
Compute a mask where smoothed residuals greater than some value
Perhaps useful for flagging contaminated pixels that aren’t in the model, such as high orders dispersed from objects that fall off of the direct image, but this hasn’t yet been extensively tested.
- Parameters
- gaussian_widthfloat
Width of the Gaussian filter used with
gaussian_filter.- thresholdfloat
Threshold, in sigma, above which to flag residuals.
- Returns
- Nothing, but pixels are masked in
self.grism.data['SCI'].
- Nothing, but pixels are masked in