MultiBeam¶
- class grizli.multifit.MultiBeam(beams, group_name=None, fcontam=0.0, psf=False, polyx=[0.3, 2.5], MW_EBV=0.0, min_mask=0.01, min_sens=0.08, sys_err=0.0, mask_resid=True, verbose=True, replace_direct=None, restore_medfilt=False, **kwargs)[source]¶
Bases:
GroupFitter
Tools for dealing with multiple
BeamCutout
instances- Parameters
- beamslist
List of
BeamCutout
objects.- group_namestr, None
Rootname to use for saved products. If None, then default to ‘group’.
- fcontamfloat
Factor to use to downweight contaminated pixels. The pixel inverse variances are scaled by the following weight factor when evaluating chi-squared of a 2D fit,
weight = np.exp(-(fcontam*np.abs(contam)*np.sqrt(ivar)))
where
contam
is the contaminating flux andivar
is the initial pixel inverse variance.- psfbool
Fit an ePSF model to the direct image to use as the morphological reference.
- polyxlist
Polynomial coefficients to use for the
compute_model
method.- MW_EBVfloat
Milky way foreground extinction.
- min_maskfloat
Minimum factor relative to the maximum pixel value of the flat f-lambda model where the 2D cutout data are considered good. Passed through to
BeamCutout
.- min_sensfloat
See
BeamCutout
.- sys_errfloat
Systematic error added in quadrature to the pixel variances:
var_total = var_initial + (beam.sci*sys_err)**2
- mask_residbool
Mask residuals above a certain threshold.
- verbosebool
Print verbose information.
- replace_directdict
Replace the direct image cutouts with a new image. Keys are
ref_image
,interp
,cutout
, andbackground_func
.- restore_medfiltbool
Restore median-filtered direct image cutouts.
- Attributes
- Asavedict
Dictionary to store the results of the fit.
Attributes Summary
Length of the beam list.
Dictionary containing number of exposures by grism.
Available PAs in each grism
Available grisms
ID of the first beam in the list.
Methods Summary
apply_trace_shift
([set_to_zero])Apply the current
beam.yoffset
to the trace of each beam.check_for_bad_PAs
([poly_order, ...])Check for bad PAs based on chi-squared values
Compute number of exposures and total exposure time for each grism.
compute_model
([id, spectrum_1d, is_cgs, ...])Compute the dispersed 2D model for an assumed input spectrum
compute_model_psf
([id, spectrum_1d, is_cgs])Compute the dispersed 2D model for an assumed input spectrum and for ePSF morphologies
drizzle_fit_lines
(fit, pline[, force_line, ...])Drizzle line maps from the individual beams.
drizzle_grisms_and_PAs
([size, fcontam, ...])Make figure showing spectra at different orients/grisms
drizzle_segmentation
([wcsobj, kernel, ...])Drizzle segmentation image from individual
MultiBeam.beams
.drizzle_segmentation_id
([id, wcsobj, ...])Drizzle segmentation image for a single ID
eval_poly_spec
(coeffs_full)Evaluate polynomial spectrum.
eval_trace_shift
(shifts, self, indices, ...)Evaluate the trace shifts.
extend
(new[, verbose])Concatenate
MultiBeam
objectsfit_at_z
([z, templates, fitter, ...])Fit the model at a given redshift.
fit_redshift
([prior, poly_order, fwhm, ...])Fit redshifts for the object.
fit_stars
([poly_order, fitter, ...])Fit stellar templates to the data.
fit_trace_shift
([split_groups, max_shift, ...])Fit a global trace shift to the beams.
flag_with_drizzled
(hdul[, sigma, update, ...])Update
MultiBeam
masks based on the blotted drizzled combined imageinit_poly_coeffs
([flat, poly_order])Initialize polynomial coefficients for the continuum fit.
load_beam_fits
(beam_list[, conf, verbose])Load a list of individual beam cutouts.
load_master_fits
(beam_file[, verbose])Load a "beams.fits" file.
Make a`~astropy.io.fits.HDUList` object with just a simple PrimaryHDU
oned_spectrum
([tfit, get_contam, ...])Compute full 1D spectrum with optional best-fit model
oned_spectrum_to_hdu
([sp, outputfile])Generate 1D spectra fits HDUList
parse_fit_outputs
(z, templates, coeffs_full, A)Parse output from
fit_at_z
.redshift_fit_twod_figure
(fit[, ...])Make figure of 2D spectrum
replace_direct_image_cutouts
([ref_image, ...])Replace "REF" extensions in a
beams.fits
filereplace_segmentation_image_cutouts
([ref_image])Replace "REF" extensions in a
beams.fits
filereshape_flat
(flat_array)Reshape a flattened array into a list of 2D arrays.
run_full_diagnostics
([pzfit, pspec2, pline, ...])Run full redshift fit and diagnostics on a
GroupFLT
object.run_individual_fits
([z, templates])Run template fits on each exposure individually to evaluate variance in line and continuum fits.
show_redshift_fit
(fit_data[, plot_flambda, ...])Show redshift fit.
write_beam_fits
([verbose])Write individual beam cutouts to FITS files.
write_master_fits
([verbose, get_hdu, strip, ...])Store all beams in a single HDU
Attributes Documentation
- N¶
Length of the beam list.
- Ngrism¶
Dictionary containing number of exposures by grism.
- PA¶
Available PAs in each grism
- grisms¶
Available grisms
- id¶
ID of the first beam in the list.
Methods Documentation
- apply_trace_shift(set_to_zero=False)[source]¶
Apply the current
beam.yoffset
to the trace of each beam.- Parameters
- set_to_zerobool
Set beam.yoffset to zero
- check_for_bad_PAs(poly_order=1, chi2_threshold=1.5, fit_background=True, reinit=True)[source]¶
Check for bad PAs based on chi-squared values
- Parameters
- poly_orderint
Polynomial order for the background fit.
- chi2_thresholdfloat
Threshold for chi-squared ratio to the minimum value.
- fit_backgroundbool
Fit a polynomial background.
- reinitbool
Reinitialize the
MultiBeam
object with thes beams that pass the threshold.
- Returns
- fit_logdict
Log of the fit results.
- keep_dictdict
Dictionary of PAs to keep.
- has_badbool
True if any PAs are flagged as bad.
- compute_model(id=None, spectrum_1d=None, is_cgs=False, apply_sensitivity=True, scale=None, reset=True)[source]¶
Compute the dispersed 2D model for an assumed input spectrum
This is a wrapper around the
grizli.model.GrismDisperser.compute_model
method.Nothing returned, but the
model
andmodelf
attributes are updated on theGrismDisperser
subcomponents of thebeams
list.- Parameters
- idstr
Object ID.
- spectrum_1darray-like
1D spectrum to use for the model.
- is_cgsbool
Spectrum is in cgs units.
- apply_sensitivitybool
Apply the sensitivity function to the model.
- scalefloat
Scale factor for the model.
- resetbool
Reset the model to zero before computing the new model.
- compute_model_psf(id=None, spectrum_1d=None, is_cgs=False)[source]¶
Compute the dispersed 2D model for an assumed input spectrum and for ePSF morphologies
This is a wrapper around the
grizli.model.GrismDisperser.compute_model_psf
method.Nothing returned, but the
model
andmodelf
attributes are updated on theGrismDisperser
subcomponents of thebeams
list.- Parameters
- idstr
Object ID.
- spectrum_1darray-like
1D spectrum to use for the model.
- is_cgsbool
Spectrum is in cgs units.
- drizzle_fit_lines(fit, pline, force_line=['Ha+NII', 'Ha', 'OIII', 'Hb', 'OII'], save_fits=True, mask_lines=True, mask_sn_limit=3, mask_4959=True, verbose=True, include_segmentation=True, get_ir_psfs=True, min_line_sn=4)[source]¶
Drizzle line maps from the individual beams.
- Parameters
- fitdict
Fit data generated by
fit_redshift
.- plinelist
List of line names to drizzle.
- force_linelist
List of line names to always drizzle.
- save_fitsbool
Save the drizzled line maps to FITS files.
- mask_linesbool
Mask pixels around other lines.
- mask_sn_limitfloat
Mask pixels around lines with S/N below this limit.
- mask_4959bool
Mask OIII-4959 if detected.
- verbosebool
Print status messages.
- include_segmentationbool
Drizzle the segmentation image.
- get_ir_psfsbool
Get the IR PSFs from the
MultiBeam.beams
.- min_line_snfloat
Minimum S/N to drizzle a line.
- Returns
- hdu_fulllist
List of
HDUList
objects with the drizzled line maps.
- drizzle_grisms_and_PAs(size=10, fcontam=0, flambda=False, scale=1, pixfrac=0.5, kernel='square', usewcs=False, tfit=None, diff=True, grism_list=['G800L', 'G102', 'G141', 'GR150C', 'GR150R', 'F090W', 'F115W', 'F150W', 'F200W', 'F356W', 'F300M', 'F335M', 'F360M', 'F410M', 'F430M', 'F460M', 'F480M', 'F444W'], mask_segmentation=True, reset_model=True, make_figure=True, fig_args={'average_only': False, 'cmap': 'viridis_r', 'mask_segmentation': True, 'scale_size': 1}, **kwargs)[source]¶
Make figure showing spectra at different orients/grisms
- Parameters
- sizeint
Image size in pixels
- fcontamfloat
Contamination parameter for drizzle.
- flambdabool
Convert to f-lambda units.
- scalefloat
Scale factor for drizzle.
- pixfracfloat
Drizzle pixfrac parameter.
- kernelstr
Drizzle kernel parameter.
- usewcsbool
Use WCS information for drizzle.
- tfitdict
Template fit results.
- diffbool
Show difference from the model.
- grism_listlist
List of grisms to include.
- mask_segmentationbool
Mask segmentation regions.
- reset_modelbool
Reset the model to a flat spectrum.
- make_figurebool
Make the figure.
- fig_argsdict
Arguments for the figure.
- Returns
- fig
Figure
Figure object.
- fig
- drizzle_segmentation(wcsobj=None, kernel='square', pixfrac=1, verbose=False)[source]¶
Drizzle segmentation image from individual
MultiBeam.beams
.- Parameters
- wcsobj: `~astropy.wcs.WCS` or `~astropy.io.fits.Header`
Output WCS.
- kernel: e.g., ‘square’, ‘point’, ‘gaussian’
Drizzle kernel, see
drizzle
.- pixfrac: float
Drizzle ‘pixfrac’, see
drizzle
.- verbose: bool
Print status messages.
- Returns
- drizzle_segmentation_id(id=None, wcsobj=None, kernel='square', pixfrac=1, verbose=True)[source]¶
Drizzle segmentation image for a single ID
- Parameters
- Returns
- eval_poly_spec(coeffs_full)[source]¶
Evaluate polynomial spectrum.
- Parameters
- coeffs_fullarray-like
Full list of polynomial coefficients.
- Returns
- xspec, yfullarray-like
Arrays of wavelength and polynomial spectrum.
- static eval_trace_shift(shifts, self, indices, poly_order, lm, verbose, fit_with_psf)[source]¶
Evaluate the trace shifts.
- Parameters
- Returns
- chi2float
Chi-squared value.
- extend(new, verbose=True)[source]¶
Concatenate
MultiBeam
objects- Parameters
- new
MultiBeam
Beam object containing new beams to add.
- verbosebool
Print summary of the change.
- new
- fit_at_z(z=0.0, templates={}, fitter='nnls', fit_background=True, poly_order=0)[source]¶
Fit the model at a given redshift.
- Parameters
- zfloat
Redshift to fit.
- templatesdict
Dictionary of template spectra to fit.
- fitterstr
Fitting method. Options are ‘nnls’ and ‘lstsq’.
- fit_backgroundbool
Fit a constant background.
- poly_orderint
Order of the polynomial to fit.
- Returns
- out_coeffsarray-like
Output coefficients of the fit.
- fit_redshift(prior=None, poly_order=1, fwhm=1200, make_figure=True, zr=None, dz=None, verbose=True, fit_background=True, fitter='nnls', delta_chi2_threshold=0.004, zoom=True, line_complexes=True, templates={}, figsize=[8, 5], fsps_templates=False)[source]¶
Fit redshifts for the object.
- Parameters
- prior: list or array-like or None
- A list or array-like object containing two elements:
prior[0]: The x-coordinates for interpolation.
prior[1]: The y-coordinates for interpolation.
Used to compute interpolated values over the
zgrid
usingnumpy.interp
when notNone
.
- poly_orderint
Polynomial order to fit.
- fwhmfloat
FWHM of the line spread function.
- make_figurebool
Generate a figure of the fit.
- zrlist or None
Redshift range to fit.
- dzlist or None
Redshift step size.
- verbosebool
Print information to the terminal.
- fit_backgroundbool
Fit a constant background.
- fitterstr
Fitting method. Options are ‘nnls’ and ‘lstsq’.
- delta_chi2_thresholdfloat
Threshold for computing the redshift grid. Not used.
- zoombool
Zoom in on the best-fit redshift.
- line_complexesbool
Include line complexes in the fit.
- templatesdict
Dictionary of template spectra to fit.
- figsizelist
Figure size for the output plot.
- fsps_templatesbool
Use FSPS templates.
- Returns
- fit_datadict
Fit data.
- fig
Figure
Figure object.
- fit_stars(poly_order=1, fitter='nnls', fit_background=True, verbose=True, make_figure=True, fwhm=0, **kwargs)[source]¶
Fit stellar templates to the data.
- Parameters
- poly_orderint
Polynomial order to fit.
- fitterstr
Fitting method. Options are ‘nnls’ and ‘lstsq’.
- fit_backgroundbool
Fit a constant background.
- verbosebool
Print information to the terminal.
- make_figurebool
Generate a figure of the fit.
- zoombool
Zoom in on the best-fit redshift. Not used.
- delta_chi2_thresholdfloat
Threshold for computing the redshift grid. Not used.
- zrfloat
Redshift to fit. Not used.
- dzfloat
Redshift step size. Not used.
- fwhmfloat
FWHM of the line spread function.
- priorNone
Not used.
- templatesdict
Dictionary of template spectra to fit.
- figsizelist
Figure size for the output plot. Not used.
- fsps_templatesbool
Not used.
- Returns
- fit_datadict
Fit data.
- fig
Figure
Figure object.
- fit_trace_shift(split_groups=True, max_shift=5, tol=0.01, verbose=True, lm=False, fit_with_psf=False, reset=False)[source]¶
Fit a global trace shift to the beams.
- Parameters
- split_groupsbool
Split the fit into groups.
- max_shiftfloat
Maximum allowed shift in pixels.
- tolfloat
Tolerance for the fit.
- verbosebool
Print status messages.
- lmbool
Use Levenberg-Marquardt optimization.
- fit_with_psfbool
Fit with the PSF model.
- resetbool
Reset the trace shifts to zero.
- Returns
- flag_with_drizzled(hdul, sigma=4, update=True, interp='nearest', verbose=True)[source]¶
Update
MultiBeam
masks based on the blotted drizzled combined image[in progress … xxx]
- Parameters
- hdul
HDUList
FITS HDU list output from
drizzle_grisms_and_PAs
or read from astack.fits
file.- sigmafloat
Residual threshold to flag.
- updatebool
Update the mask.
- interpstr
Interpolation method for
ablot
.- verbosebool
Print status messages.
- hdul
- Returns
- Updates the individual
fit_mask
attributes of the individual beams - if
update==True
.
- Updates the individual
- init_poly_coeffs(flat=None, poly_order=1)[source]¶
Initialize polynomial coefficients for the continuum fit.
- Parameters
- flatarray-like
Flattened array of the flat f-lambda model.
- poly_orderint
Order of the polynomial to fit.
- Returns
- None (updates
self.A_poly
,self.n_poly
,self.x_poly
).
- None (updates
- load_beam_fits(beam_list, conf=None, verbose=True)[source]¶
Load a list of individual beam cutouts.
- Parameters
- beam_listlist
List of filenames to load.
- conf
GroupFitter
Configuration object to use for the beam cutouts.
- verbosebool
Print the name of the file being loaded.
- load_master_fits(beam_file, verbose=True)[source]¶
Load a “beams.fits” file.
- Parameters
- beam_filestr
Filename of the “beams.fits” file.
- verbosebool
Print the name of the file being loaded.
- make_simple_hdulist()[source]¶
Make a`~astropy.io.fits.HDUList` object with just a simple PrimaryHDU
- oned_spectrum(tfit=None, get_contam=True, get_background=False, masked_model=None, **kwargs)[source]¶
Compute full 1D spectrum with optional best-fit model
- oned_spectrum_to_hdu(sp=None, outputfile=None, **kwargs)[source]¶
Generate 1D spectra fits HDUList
- Parameters
- spoptional, dict
Output of
oned_spectrum
. If None, then run that function with**kwargs
.- outputfileNone, str
If a string supplied, then write the
HDUList
to a file.
- Returns
- hdul
HDUList
FITS version of the 1D spectrum tables.
- hdul
- parse_fit_outputs(z, templates, coeffs_full, A)[source]¶
Parse output from
fit_at_z
.- Parameters
- zfloat
Redshift at which to evaluate the fits.
- templateslist of
SpectrumTemplate
objects Generated with, e.g.,
load_templates
.- coeffs_full
ndarray
Template fit coefficients
- A
ndarray
Matrix generated for fits and used for computing model 2D spectra:
>>> model_flat = np.dot(coeffs_full, A) >>> # mb = MultiBeam(...) >>> all_models = mb.reshape_flat(model_flat) >>> m0 = all_models[0] # model for mb.beams[0]
- Returns
- line_fluxdict
Line fluxes and uncertainties, in cgs units (erg/s/cm2)
- covar
ndarray
Covariance matrix for the fit coefficients
- cont1d, line1d, model1d
SpectrumTemplate
Best-fit continuum, line, and full (continuum + line) templates
- model_continuum
ndarray
Flat array of the best fit 2D continuum
- redshift_fit_twod_figure(fit, spatial_scale=1, dlam=46.0, NY=10, figsize=[8, 3.5], **kwargs)[source]¶
Make figure of 2D spectrum
- Parameters
- fitdict
Fit data generated by
fit_redshift
.- spatial_scalefloat
Scale factor for the spatial axis.
- dlamfloat
Pixel size in Angstroms.
- NYint
Number of pixels in the spatial direction.
- figsizelist
Figure size for the output plot.
- Returns
- replace_direct_image_cutouts(ref_image='gdn-100mas-f160w_drz_sci.fits', ext=0, interp='poly5', cutout=200, background_func=<function median>, thumb_labels=None)[source]¶
Replace “REF” extensions in a
beams.fits
file- Parameters
- ref_imagestr or
HDUList
Filename or preloaded FITS file.
- extint
Extension of the
ref_image
to use.- interpstr
Interpolation function to use for
do_blot
.- cutoutint
Make a slice of the
ref_image
with size [-cutout,+cutout] around the center position of the desired object before passing toblot
.- background_funcfunction, None
If not
None
, compute local background with value frombackground_func(ref_image[cutout])
.- thumb_labelslist
List of thumbnail labels to replace.
- ref_imagestr or
- Returns
- beams_image
HDUList
Image object with the “REF” extensions filled with the new blotted image cutouts.
- beams_image
- replace_segmentation_image_cutouts(ref_image='gdn-100mas-f160w_seg.fits')[source]¶
Replace “REF” extensions in a
beams.fits
file
- reshape_flat(flat_array)[source]¶
Reshape a flattened array into a list of 2D arrays.
- Parameters
- flat_arrayarray-like
Flattened array to reshape.
- Returns
- outlist
List of 2D arrays.
- run_full_diagnostics(pzfit={}, pspec2={}, pline={}, force_line=['Ha+NII', 'Ha', 'OIII', 'Hb', 'OII'], GroupFLT=None, prior=None, zoom=True, verbose=True)[source]¶
Run full redshift fit and diagnostics on a
GroupFLT
object.- Parameters
- pzfitdict
Parameters for the redshift fit.
- pspec2dict
Parameters for the 2D spectrum drizzle.
- plinedict
Parameters for the line drizzle.
- force_linelist
List of line names to always drizzle.
- GroupFLT
GroupFLT
Group object with the full set of objects to fit simultaneously.
- prior
Table
Prior catalog for the redshift fit.
- zoombool
Zoom in on the object.
- verbosebool
Print status messages.
- Returns
- run_individual_fits(z=0, templates={})[source]¶
Run template fits on each exposure individually to evaluate variance in line and continuum fits.
- Parameters
- zfloat
Redshift at which to evaluate the fit
- templateslist of
SpectrumTemplate
objects Generated with, e.g.,
load_templates
.
- Returns
- line_flux, line_errdict
Dictionaries with the measured line fluxes and uncertainties for each exposure fit.
- coeffs_list
ndarray
[Nbeam x Ntemplate] Raw fit coefficients
- chi2_list, DoF_list
ndarray
[Nbeam] Chi-squared and effective degrees of freedom for each separate fit
- show_redshift_fit(fit_data, plot_flambda=True, figsize=[8, 5])[source]¶
Show redshift fit.
- Parameters
- fit_datadict
Fit data generated by
fit_redshift
.- plot_flambdabool
Plot in f-lambda units.
- figsizelist
Figure size for the output plot.
- Returns
- fig
Figure
Figure object.
- fig
- write_beam_fits(verbose=True)[source]¶
Write individual beam cutouts to FITS files.
- Parameters
- verbosebool
Print the name of the output file.
- Returns
- outfileslist
List of output filenames.