run_all

grizli.fitting.run_all(id, t0=None, t1=None, fwhm=1200, zr=[0.65, 1.6], dz=[0.004, 0.0002], fitter=['nnls', 'bounded'], group_name='grism', fit_stacks=True, only_stacks=False, prior=None, fcontam=0.2, pline={'kernel': 'point', 'pixfrac': 0.2, 'pixscale': 0.1, 'size': 8, 'wcs': None}, min_line_sn=4, mask_sn_limit=inf, fit_only_beams=False, fit_beams=True, root='*', fit_trace_shift=False, phot=None, use_phot_obj=True, phot_obj=None, verbose=True, scale_photometry=False, show_beams=True, scale_on_stacked_1d=True, use_cached_templates=True, loglam_1d=True, overlap_threshold=5, MW_EBV=0.0, sys_err=0.03, huber_delta=4, get_student_logpdf=False, get_dict=False, bad_pa_threshold=1.6, units1d='flam', redshift_only=False, line_size=1.6, dscale=0.25, scale_linemap=1, use_psf=False, get_line_width=False, sed_args={'bin': 1, 'xlim': [0.3, 9]}, get_ir_psfs=True, min_mask=0.01, min_sens=0.02, mask_resid=True, save_stack=True, full_line_list=['Lya', 'OII', 'Hb', 'OIII', 'Ha', 'Ha+NII', 'SII', 'SIII'], get_line_deviations=True, bounded_kwargs={'method': 'bvls', 'tol': 1e-08, 'verbose': 0}, write_fits_files=True, save_figures=True, fig_type='png', diff2d=True, **kwargs)[source]

Run the full template-fitting procedure

  1. Load MultiBeam and stack files

  2. … tbd

Parameters
idint

Object ID in the internal catalogs. This is generally an int, but in principle could be a str or something else.

t0dict

Dictionary of SpectrumTemplate objects used for the redshift fits. Generally these will have fixed line ratios to avoid unphysical line degeneracies (e.g., very strong [SII] without H-alpha).

If None, then the templates are generated with

>>> t0 = grizli.utils.load_templates(line_complexes=True, fsps_templates=True, fwhm=fwhm)
t1dict

Dictionary of SpectrumTemplate objects used for the final fit at the best-fit redshift. Generally these will be separate continuum and individual line templates so that the line fluxes are determined freely (which are then also needed if you want to make the drizzled narrowband emission line maps).

If None, then the templates are generated with

>>> t1 = grizli.utils.load_templates(line_complexes=False, fsps_templates=True, fwhm=fwhm)

Note

As of 66a3ec5 all templates can be eazy.templates.Template objects.

fwhmfloat

Line FWHM passed to load_templates if t0 or t1 not specified.

zr[float, float], [float], or 0

Redshift range to fit.

  • [z1, z2] - fit on a logarithmic grid between z1 and z2 with steps specified in dz

  • [zfix] - fit templates at a specified value

  • 0 - fit stellar templates only

dz[float, float]

Logarithmic step size (1+z) of redshift grid. See log_zgrid.

fitter[str, str]

Least squares optimization method (‘nnls’,’lstsq’,’bounded’). The first option is used for the redshift fit with the t0 templates and the second is used for the final fit with the t1 templates.

  • nnls: Generally SPS continuum templates should be fit with nnls to enforce physical template combinations.

  • bounded: Enforces non-negative continuum templates but allows line templates (with a name starting with line [space]) to be negative. The bounded fits are controlled with bounded_kwargs and the flux limits set in the global parameter grizli.fitting.LINE_BOUNDS.

  • lstsq: Fit with regular least squares, e.g., for PCA templates that can have negative coefficients (e.g., load_sdss_pca_templates).

bounded_kwargsdict

Keywords passed to scipy.optimize.lsq_linear for ‘bounded’ fits.

group_namestr

Passed to MultiBeam on initialization

rootstr

Basename MultiBeam FITS filenames to search for, e.g., to concatenate separate G141 and G102 files of a single object:

>>> mb_files = glob.glob(f'{root}_{id:05d}.beams.fits')
fit_stacksbool

Fit redshifts on the stacked spectra, which can be much faster than for the separate “beams” fits, but where the model generation isn’t as robust. This is generally deprecated, but should still run.

only_stacksbool

Only fit the stacks.

priorNone, (array, array)

Redshift prior (z, pz) passed to xfit_redshift.

fcontam, min_mask, min_sens, mask_residfloat, float, float, bool

Contamination weighting passed to MultiBeam

plinedict

Parameters for drizzled line maps.

min_line_snfloat

If finite, then pass to drizzle_fit_lines to determine which line maps to create.

mask_sn_limitfloat

SN limit to pass to drizzle_fit_lines

fit_only_beamsbool

If True, only fit with MultiBeam objects.

fit_beamsbool

Fit with MultiBeam objects.

fit_trace_shiftbool

Fit for shifts of the traces fo each group oof beams.

photNone, dict

Photometry dictionary passed to set_photometry

use_phot_objbool

Use phot_obj if it is available.

phot_objNone, EazyPhot

Catalog object for automatically generating phot dictionaries

verbosebool

Some control over the runtime verbosity

scale_photometrybool

If photometry is available, try to normalize the spectra and photometry.

show_beams, scale_on_stacked_1d, loglam_1dbool, bool

Passed to xmake_fit_plot for the final redshift fit plot.

use_cached_templatesbool

Passed to xfit_at_z

overlap_thresholdfloat

Parameter for StackFitter when fitting on stacks.

MW_EBVfloat

Galactic extinction E(B-V) (mag)

sys_errfloat

Systematic error used for the spectra and photometry, multiplied to the flux densities and added in quadrature to the nominal uncertainties.

huber_deltafloat

Passed to xfit_at_z for using a Huber loss function.

get_student_logpdfbool

Use Student-t likelihood on redshift_fit

get_dictbool

Don’t actually run anything, just return a dictionary with all of the keyword parameters passed to the function

bad_pa_thresholdfloat

Threshold for identifying bad PAs when using StackFitter objects (not beams)

units1dstr

Not used

redshift_onlybool

Just run the redshift fit, don’t drizzle the line maps

line_sizefloat

Cutout size in arcsec of the line map figures.

use_psfbool

Initialize the MultiBeam objects with psf=True to fit the morphology using the EffectivePSF models.

get_line_widthbool

Try to fit for emission line velocity widths (developmental)

sed_argsdict

Keyword arguments passed to full_sed_plot when photometry + spectra are available

get_ir_psfsbool

Include PSF extensions in the drizzled line maps derived from the EffectivePSF models.

save_stackbool

Generate a stack.fits file from the beams fit

full_line_listlist

Line list passed to show_drizzled_lines to determine which lines are always included in the drizzled line maps.

get_line_deviationsbool

Check plausibility of fit coefficients with check_tfit_coeffs

write_fits_filesbool

Save ‘full.fits’ and ‘stack.fits’ files

save_figures, fig_typebool, str

Save diagnostic figure files with extension fig_type

Returns
mbMultiBeam

The beams object used for the redshift / template fits

stStackFitter

The stacked spectrum object generated from the ‘beams’

fitastropy.table.Table

Table with the fit results

tfitdict

Various parameters of the template fit at the final redshift

line_hduHDUList

Drizzled line maps