API

fastspecfit

Tools for fast stellar continuum, emission-line, and broadband photometric modeling.

fastspecfit.continuum

Methods and tools for continuum-fitting.

class fastspecfit.continuum.ContinuumTools(ignore_photometry=False, fphoto=None, emlinesfile=None)[source]

Tools for dealing with stellar continua.

Parameters:
  • templates (str, optional) – Full path to the templates used for continuum-fitting.

  • templateversion (str, optional, defaults to v1.0) – Version of the templates.

  • mapdir (str, optional) – Full path to the Milky Way dust maps.

  • mintemplatewave (float, optional, defaults to None) – Minimum template wavelength to read into memory. If None, the minimum available wavelength is used (around 100 Angstrom).

  • maxtemplatewave (float, optional, defaults to 6e4) – Maximum template wavelength to read into memory.

  • note:: (..) – Need to document all the attributes.

static build_linemask(wave, flux, ivar, redshift=0.0, nsig=7.0, linetable=None, emlinesfile=None, log=None, verbose=False)[source]

Generate a mask which identifies pixels impacted by emission lines.

Parameters:
  • wave (numpy.ndarray [npix]) – Observed-frame wavelength array.

  • flux (numpy.ndarray [npix]) – Spectrum corresponding to wave.

  • ivar (numpy.ndarray [npix]) – Inverse variance spectrum corresponding to flux.

  • redshift (float, optional, defaults to 0.0) – Object redshift.

  • nsig (float, optional, defaults to 5.0) – Mask pixels which are within +/-nsig-sigma of each emission line, where sigma is the line-width in km/s.

Returns:

Smooth continuum spectrum which can be subtracted from flux in order to create a pure emission-line spectrum.

Return type:

smooth numpy.ndarray [npix]

dict with the following keys:

linemask : list linemask_all : list linename : list linepix : list contpix : list

Notes

Code exists to generate some QA but the code is not exposed.

static call_nnls(modelflux, flux, ivar, xparam=None, debug=False, interpolate_coeff=False, xlabel=None, png=None, log=None)[source]

Wrapper on nnls.

Works with both spectroscopic and photometric input and with both 2D and 3D model spectra.

To be documented.

interpolate_coeff - return the interpolated coefficients when exploring

an array or grid of xparam

continuum_fluxes(data, templatewave, continuum, log=None)[source]

Compute rest-frame luminosities and observed-frame continuum fluxes.

kcorr_and_absmag(data, templatewave, continuum, snrmin=2.0, log=None)[source]

Compute K-corrections and absolute magnitudes.

static smooth_and_resample(templateflux, templatewave, specwave=None, specres=None)[source]

Given a single template, apply the resolution matrix and resample in wavelength.

Parameters:
  • templateflux (numpy.ndarray [npix]) – Input (model) spectrum.

  • templatewave (numpy.ndarray [npix]) – Wavelength array corresponding to templateflux.

  • specwave (numpy.ndarray [noutpix], optional, defaults to None) – Desired output wavelength array, usually that of the object being fitted.

  • specres (desispec.resolution.Resolution, optional, defaults to None) – Resolution matrix.

Returns:

Smoothed and resampled flux at the new resolution and wavelength sampling.

Return type:

numpy.ndarray [noutpix]

Notes

This function stands by itself rather than being in a class because we call it with multiprocessing, below.

templates2data(_templateflux, _templatewave, redshift=0.0, dluminosity=None, vdisp=None, cameras=['b', 'r', 'z'], specwave=None, specres=None, specmask=None, coeff=None, photsys=None, synthphot=True, stack_cameras=False, debug=False, log=None)[source]

Work-horse routine to turn input templates into spectra that can be compared to real data.

Redshift, apply the resolution matrix, and resample in wavelength.

Parameters:
  • redshift

  • specwave

  • specres

  • south

  • photometry? (synthphot - synthesize)

Return type:

Vector or 3-element list of [npix, nmodel] spectra.

Notes

This method does none or more of the following: - redshifting - wavelength resampling - apply velocity dispersion broadening - apply the resolution matrix - synthesize photometry

It also naturally handles templates which have been precomputed on a grid of velocity dispersion (and therefore have an additional dimension). However, if the input grid is 3D, it is reshaped to be 2D but then it isn’t reshaped back because of the way the photometry table is organized (bug or feature?).

fastspecfit.continuum._convolve_vdisp(templateflux, vdisp, pixsize_kms)[source]

Convolve by the velocity dispersion.

Parameters:
  • templateflux

  • vdisp

Notes

fastspecfit.continuum._estimate_linesigmas(wave, flux, ivar, redshift=0.0, png=None, refit=True, log=None, verbose=False)[source]

Estimate the velocity width from potentially strong, isolated lines.

fastspecfit.continuum._smooth_continuum(wave, flux, ivar, redshift, camerapix=None, medbin=175, smooth_window=75, smooth_step=25, maskkms_uv=3000.0, maskkms_balmer=1000.0, maskkms_narrow=200.0, linetable=None, emlinesfile=None, linemask=None, png=None, log=None, verbose=False)[source]

Build a smooth, nonparametric continuum spectrum.

Parameters:
  • wave (numpy.ndarray [npix]) – Observed-frame wavelength array.

  • flux (numpy.ndarray [npix]) – Spectrum corresponding to wave.

  • ivar (numpy.ndarray [npix]) – Inverse variance spectrum corresponding to flux.

  • redshift (float) – Object redshift.

  • medbin (int, optional, defaults to 150 pixels) – Width of the median-smoothing kernel in pixels; a magic number.

  • smooth_window (int, optional, defaults to 50 pixels) – Width of the sliding window used to compute the iteratively clipped statistics (mean, median, sigma); a magic number. Note: the nominal extraction width (0.8 A) and observed-frame wavelength range (3600-9800 A) corresponds to pixels that are 66-24 km/s. So smooth_window of 50 corresponds to 3300-1200 km/s, which is conservative for all but the broadest lines. A magic number.

  • smooth_step (int, optional, defaults to 10 pixels) – Width of the step size when computing smoothed statistics; a magic number.

  • maskkms_uv (float, optional, defaults to 3000 km/s) – Masking width for UV emission lines. Pixels within +/-3*maskkms_uv are masked before median-smoothing.

  • maskkms_balmer (float, optional, defaults to 3000 km/s) – Like maskkms_uv but for Balmer lines.

  • maskkms_narrow (float, optional, defaults to 300 km/s) – Like maskkms_uv but for narrow, forbidden lines.

  • linemask (numpy.ndarray of type bool, optional, defaults to None) – Boolean mask with the same number of pixels as wave where True means a pixel is (possibly) affected by an emission line (specifically a strong line which likely cannot be median-smoothed).

  • png (str, optional, defaults to None) – Generate a simple QA plot and write it out to this filename.

Returns:

  • smooth numpy.ndarray [npix] – Smooth continuum spectrum which can be subtracted from flux in order to create a pure emission-line spectrum.

  • smoothsigma numpy.ndarray [npix] – Smooth one-sigma uncertainty spectrum.

fastspecfit.continuum.continuum_specfit(data, result, templatecache, fphoto=None, emlinesfile=None, constrain_age=False, no_smooth_continuum=False, ignore_photometry=False, fastphot=False, log=None, debug_plots=False, verbose=False)[source]

Fit the non-negative stellar continuum of a single spectrum.

Parameters:

data (dict) – Dictionary of input spectroscopy (plus ancillary data) populated by fastspecfit.io.DESISpectra.read_and_unpack().

Returns:

Table with all the continuum-fitting results with columns documented in init_output().

Return type:

astropy.table.Table

Notes

  • Consider using cross-correlation to update the redrock redshift.

  • We solve for velocity dispersion if …

fastspecfit.continuum.restframe_photometry(redshift, zmodelflux, zmodelwave, maggies, ivarmaggies, filters_in, absmag_filters, band_shift=None, snrmin=2.0, dmod=None, cosmo=None, log=None)[source]

Compute K-corrections and rest-frame photometry for a single object.

Parameters:
  • redshift (float) – Galaxy or QSO redshift.

  • zmodelwave (numpy.ndarray) – Observed-frame (redshifted) model wavelength array.

  • zmodelflux (numpy.ndarray) – Observed-frame model spectrum.

  • maggies (numpy.ndarray) – Input photometric fluxes in the filters_in bandpasses.

  • ivarmaggies (numpy.ndarray) – Inverse variance photometry corresponding to maggies.

  • filters_in (speclite.filters.FilterSequence) – Input filter curves.

  • absmag_filters (speclite.filters.FilterSequence) – Filter curves corresponding to desired bandpasses.

  • band_shift (numpy.ndarray or None) – Band-shift each bandpass in absmag_filters by this amount.

  • snrmin (float, defaults to 2.) – Minimum signal-to-noise ratio in the input photometry (maggies) in order for that bandpass to be used to compute a K-correction.

  • dmod (float or None) – Distance modulus corresponding to redshift. Not needed if cosmo is provided.

  • cosmo (fastspecfit.util.TabulatedDESI or None) – Cosmological model class needed to compute the distance modulus.

  • log (desiutil.log) – Logging object.

Returns:

  • kcorr (numpy.ndarray) – K-corrections for each bandpass in absmag_filters.

  • absmag (numpy.ndarray) – Absolute magnitudes, band-shifted according to band_shift (if provided) for each bandpass in absmag_filters.

  • ivarabsmag (numpy.ndarray) – Inverse variance corresponding to absmag.

  • synth_absmag (numpy.ndarray) – Like absmag, but entirely based on synthesized photometry.

  • synth_maggies_in (numpy.ndarray) – Synthesized input photometry (should closely match maggies if the model fit is good).

Notes

By default, the K-correction is computed by finding the observed-frame bandpass closest in wavelength (and with a minimum signal-to-noise ratio) to the desired band-shifted absolute magnitude bandpass. In other words, by default we endeavor to minimize the K-correction. The inverse variance, ivarabsmag, is derived from the inverse variance of the K-corrected photometry. If no bandpass is available then ivarabsmag is set to zero and absmag is derived from the synthesized rest-frame photometry.

fastspecfit.emlines

Methods and tools for fitting emission lines.

fastspecfit.emlines._objective_function(free_parameters, emlinewave, emlineflux, weights, redshift, dlog10wave, resolution_matrix, camerapix, parameters, Ifree, Itied, tiedtoparam, tiedfactor, doubletindx, doubletpair, linewaves)[source]

The parameters array should only contain free (not tied or fixed) parameters.

fastspecfit.emlines.build_emline_model(dlog10wave, redshift, lineamps, linevshifts, linesigmas, linewaves, emlinewave, resolution_matrix, camerapix=None, debug=False)[source]

Given parameters, build the model emission-line spectrum.

ToDo: can this be optimized using numba?

fastspecfit.emlines.emline_specfit(data, templatecache, result, continuummodel, smooth_continuum, minsnr_balmer_broad=3.0, fphoto=None, emlinesfile=None, synthphot=True, broadlinefit=True, percamera_models=False, log=None, verbose=False)[source]

Perform the fit minimization / chi2 minimization.

Parameters:
  • data

  • continuummodel

  • smooth_continuum

  • synthphot

  • verbose

  • broadlinefit

Returns:

  • results

  • modelflux

fastspecfit.emlines.read_emlines(emlinesfile=None)[source]

Read the set of emission lines of interest.

fastspecfit.fastspecfit

See sandbox/running-fastspecfit for examples.

fastspecfit.fastspecfit._assign_units_to_columns(fastfit, metadata, Spec, templates, fastphot, stackfit, log)[source]

Assign astropy units to output tables.

fastspecfit.fastspecfit._fastspec_one(args)[source]

Multiprocessing wrapper.

fastspecfit.fastspecfit.fastphot(args=None, comm=None)[source]

Main fastphot script.

This script is the engine to model the broadband photometry of one or more DESI objects. It initializes the ContinuumFit class, reads the data, fits each object (with the option of fitting in parallel), and writes out the results as a multi-extension binary FITS table.

Parameters:
  • args (argparse.Namespace or None) – Required and optional arguments parsed via inputs to the command line.

  • comm (mpi4py.MPI.MPI.COMM_WORLD or None) – Intracommunicator used with MPI parallelism.

fastspecfit.fastspecfit.fastspec(fastphot=False, stackfit=False, args=None, comm=None, verbose=False)[source]

Main fastspec script.

This script is the engine to model one or more DESI spectra. It initializes the FastFit class, reads the data, fits each spectrum (with the option of fitting in parallel), and writes out the results as a multi-extension binary FITS table.

Parameters:
  • args (argparse.Namespace or None) – Required and optional arguments parsed via inputs to the command line.

  • comm (mpi4py.MPI.MPI.COMM_WORLD or None) – Intracommunicator used with MPI parallelism.

fastspecfit.fastspecfit.fastspec_one(iobj, data, out, meta, fphoto, templates, log=None, emlinesfile=None, broadlinefit=True, fastphot=False, constrain_age=False, no_smooth_continuum=False, ignore_photometry=False, percamera_models=False, debug_plots=False, minsnr_balmer_broad=3.0)[source]

Multiprocessing wrapper to run fastspec() on a single object.

fastspecfit.fastspecfit.parse(options=None, log=None)[source]

Parse input arguments to fastspec and fastphot scripts.

fastspecfit.fastspecfit.stackfit(args=None, comm=None)[source]

Wrapper script to fit (model) generic stacked spectra.

Parameters:
  • args (argparse.Namespace or None) – Required and optional arguments parsed via inputs to the command line.

  • comm (mpi4py.MPI.MPI.COMM_WORLD or None) – Intracommunicator used with MPI parallelism.

fastspecfit.html

Code to support building HTML visualizations.

fastspecfit.html._build_htmlpage_one(args)[source]

Wrapper function for the multiprocessing.

fastspecfit.html.build_htmlhome(htmldir, fastfit, metadata, tileinfo, coadd_type='cumulative', specprod='denali', htmlhome='index.html', mp=1)[source]

Build the home (index.html) page.

fastspecfit.html.build_htmlpage_one(htmlhome, survey, targetclass, zcolumn, meta, httpfile, pngfile, targiddatafile, nexttargetid, prevtargetid, nexttargidhtmlfile, prevtargidhtmlfile)[source]

Build the web page for a single object.

fastspecfit.io

Tools for reading DESI spectra and reading and writing fastspecfit files.

class fastspecfit.io.DESISpectra(stackfit=False, redux_dir=None, fiberassign_dir=None, fphotodir=None, fphotofile=None, mapdir=None)[source]
_gather_photometry(specprod=None, alltiles=None)[source]

Gather photometry. Unfortunately some of the bandpass information here will be repeated (and has to be consistent with) continuum.Fiters.

read_and_unpack(fastphot=False, synthphot=True, ignore_photometry=False, verbose=False, mp=1)[source]

Read and unpack selected spectra or broadband photometry.

Parameters:
  • fastphot (bool) – Read and unpack the broadband photometry; otherwise, handle the DESI three-camera spectroscopy. Optional; defaults to False.

  • synthphot (bool) – Synthesize photometry from the coadded optical spectrum. Optional; defaults to True.

  • remember_coadd (bool) – Add the coadded spectrum to the returned dictionary. Optional; defaults to False (in order to reduce memory usage).

Returns:

targetidnumpy.int64

DESI target ID.

zredrocknumpy.float64

Redrock redshift.

cameraslist

List of camera names present for this spectrum.

wavelist

Three-element list of numpy.ndarray wavelength vectors, one for each camera.

fluxlist

Three-element list of numpy.ndarray flux spectra, one for each camera and corrected for Milky Way extinction.

ivarlist

Three-element list of numpy.ndarray inverse variance spectra, one for each camera.

reslist

Three-element list of desispec.resolution.Resolution objects, one for each camera.

snrnumpy.ndarray

Median per-pixel signal-to-noise ratio in the grz cameras.

linemasklist

Three-element list of numpy.ndarray boolean emission-line masks, one for each camera. This mask is used during continuum-fitting.

linenamelist

Three-element list of emission line names which might be present in each of the three DESI cameras.

linepixlist

Three-element list of pixel indices, one per camera, which were identified in FFit.build_linemask to belong to emission lines.

contpixlist

Three-element list of pixel indices, one per camera, which were identified in FFit.build_linemask to not be “contaminated” by emission lines.

coadd_wavenumpy.ndarray

Coadded wavelength vector with all three cameras combined.

coadd_fluxnumpy.ndarray

Flux corresponding to coadd_wave.

coadd_ivarnumpy.ndarray

Inverse variance corresponding to coadd_flux.

photsysstr

Photometric system.

photastropy.table.Table

Total photometry in grzW1W2, corrected for Milky Way extinction.

fiberphotastropy.table.Table

Fiber photometry in grzW1W2, corrected for Milky Way extinction.

fibertotphotastropy.table.Table

Fibertot photometry in grzW1W2, corrected for Milky Way extinction.

synthphotastropy.table.Table

Photometry in grz synthesized from the Galactic extinction-corrected coadded spectra (with a mild extrapolation of the data blueward and redward to accommodate the g-band and z-band filter curves, respectively.

Return type:

List of dictionaries (dict, one per object) the following keys

read_stacked(stackfiles, firsttarget=0, ntargets=None, stackids=None, synthphot=True, ignore_photometry=False, mp=1)[source]

Read one or more stacked spectra.

Parameters:
  • stackfiles (str or array) – Full path to one or more input stacked-spectra file(s).

  • stackids (int or array or None) – Restrict the sample to the set of stackids in this list. If None, fit everything.

  • firsttarget (int) – Integer offset of the first object to consider in each file. Useful for debugging and testing. Defaults to 0.

  • ntargets (int or None) – Number of objects to analyze in each file. Useful for debugging and testing. If None, select all targets which satisfy the selection criteria.

  • synthphot (bool) – Synthesize photometry from the coadded optical spectrum. Optional; defaults to True.

Returns:

targetidnumpy.int64

DESI target ID.

zredrocknumpy.float64

Redrock redshift.

cameraslist

List of camera names present for this spectrum.

wavelist

Three-element list of numpy.ndarray wavelength vectors, one for each camera.

fluxlist

Three-element list of numpy.ndarray flux spectra, one for each camera and corrected for Milky Way extinction.

ivarlist

Three-element list of numpy.ndarray inverse variance spectra, one for each camera.

reslist

Three-element list of desispec.resolution.Resolution objects, one for each camera.

snrnumpy.ndarray

Median per-pixel signal-to-noise ratio in the grz cameras.

linemasklist

Three-element list of numpy.ndarray boolean emission-line masks, one for each camera. This mask is used during continuum-fitting.

linenamelist

Three-element list of emission line names which might be present in each of the three DESI cameras.

linepixlist

Three-element list of pixel indices, one per camera, which were identified in FFit.build_linemask to belong to emission lines.

contpixlist

Three-element list of pixel indices, one per camera, which were identified in FFit.build_linemask to not be “contaminated” by emission lines.

coadd_wavenumpy.ndarray

Coadded wavelength vector with all three cameras combined.

coadd_fluxnumpy.ndarray

Flux corresponding to coadd_wave.

coadd_ivarnumpy.ndarray

Inverse variance corresponding to coadd_flux.

photsysstr

Photometric system.

photastropy.table.Table

Total photometry in grzW1W2, corrected for Milky Way extinction.

fiberphotastropy.table.Table

Fiber photometry in grzW1W2, corrected for Milky Way extinction.

fibertotphotastropy.table.Table

Fibertot photometry in grzW1W2, corrected for Milky Way extinction.

synthphotastropy.table.Table

Photometry in grz synthesized from the Galactic extinction-corrected coadded spectra (with a mild extrapolation of the data blueward and redward to accommodate the g-band and z-band filter curves, respectively.

Return type:

List of dictionaries (dict, one per object) the following keys

static resolve(targets)[source]

Resolve which targets are primary in imaging overlap regions.

Parameters:

targets (ndarray) – Rec array of targets. Must have columns “RA” and “DEC” and either “RELEASE” or “PHOTSYS” or “TARGETID”.

Returns:

The original target list trimmed to only objects from the “northern” photometry in the northern imaging area and objects from “southern” photometry in the southern imaging area.

Return type:

ndarray

select(redrockfiles, zmin=None, zmax=None, zwarnmax=None, targetids=None, firsttarget=0, ntargets=None, input_redshifts=None, specprod_dir=None, use_quasarnet=True, redrockfile_prefix='redrock-', specfile_prefix='coadd-', qnfile_prefix='qso_qn-')[source]

Select targets for fitting and gather the necessary spectroscopic metadata.

Parameters:
  • redrockfiles (str or array) – Full path to one or more input Redrock file(s).

  • zmin (float) – Minimum redshift of observed targets to select. Defaults to 0.001. Note that any value less than or equal to zero will raise an exception because a positive redshift is needed to compute the distance modulus when modeling the broadband photometry.

  • zmax (float or None) – Maximum redshift of observed targets to select. None is equivalent to not having an upper redshift limit.

  • zwarnmax (int or None) – Maximum Redrock zwarn value for selected targets. None is equivalent to not having any cut on zwarn.

  • targetids (int or array or None) – Restrict the sample to the set of targetids in this list. If None, select all targets which satisfy the selection criteria.

  • firsttarget (int) – Integer offset of the first object to consider in each file. Useful for debugging and testing. Defaults to 0.

  • ntargets (int or None) – Number of objects to analyze in each file. Useful for debugging and testing. If None, select all targets which satisfy the selection criteria.

  • input_redshifts (float or array or None) – Input redshifts to use for each input targetids input If None, use the nominal Redrock (or QuasarNet) redshifts.

  • use_quasarnet (bool) – Use QuasarNet to improve QSO redshifts, if the afterburner file is present. Defaults to True.

  • redrockfile_prefix (str) – Prefix of the redrockfiles. Defaults to redrock-.

  • specfile_prefix (str) – Prefix of the spectroscopic coadds corresponding to the input Redrock file(s). Defaults to coadd-.

  • qnfile_prefix (str) – Prefix of the QuasarNet afterburner file. Defaults to qso_qn-.

coadd_type

Type of coadded spectra (healpix, cumulative, pernight, or perexp).

Type:

str

meta

Array of tables (one per input redrockfile) with the metadata needed to fit the data and to write to the output file(s).

Type:

list of astropy.table.Table

redrockfiles

Input Redrock file names.

Type:

str array

specfiles

Spectroscopic file names corresponding to each each input Redrock file.

Type:

str array

specprod

Spectroscopic production name for the input Redrock file.

Type:

str

Notes

We assume that specprod is the same for all input Redrock files, although we don’t explicitly do this check. Specifically, we only read the header of the first file.

fastspecfit.io._unpack_one_spectrum(args)[source]

Multiprocessing wrapper.

fastspecfit.io._unpack_one_stacked_spectrum(args)[source]

Multiprocessing wrapper.

fastspecfit.io.cache_templates(templates=None, templateversion='1.3.0', imf='chabrier', mintemplatewave=None, maxtemplatewave=400000.0, vdisp_nominal=125.0, read_linefluxes=False, fastphot=False, log=None)[source]

“Read the templates into a dictionary.

fastspecfit.io.get_qa_filename(metadata, coadd_type, outprefix=None, outdir=None, fastphot=False, log=None)[source]

Build the QA filename.

fastspecfit.io.get_templates_filename(templateversion='1.3.0', imf='chabrier')[source]

Get the templates filename.

fastspecfit.io.init_fastspec_output(input_meta, specprod, fphoto=None, templates=None, ncoeff=None, data=None, log=None, fastphot=False, emlinesfile=None, stackfit=False)[source]

Initialize the fastspecfit output data and metadata table.

Parameters:
  • tile (str) – Tile number.

  • night (str) – Night on which tile was observed.

  • redrock (astropy.table.Table) – Redrock redshift table (row-aligned to fibermap).

  • fibermap (astropy.table.Table) – Fiber map (row-aligned to redrock).

Notes

Must provide templates or ncoeff.

fastspecfit.io.one_desi_spectrum(survey, program, healpix, targetid, specprod='fuji', outdir='.', overwrite=False)[source]

Utility function to write a single DESI spectrum (e.g., for paper figures or unit tests).

fastspecfit.io.read_fastspecfit(fastfitfile, rows=None, columns=None, read_models=False)[source]

Read the fitting results.

fastspecfit.io.select(fastfit, metadata, coadd_type, healpixels=None, tiles=None, nights=None, return_index=False)[source]

Optionally trim to a particular healpix or tile and/or night.

fastspecfit.io.unpack_one_spectrum(iobj, specdata, meta, ebv, fphoto, fastphot, synthphot, ignore_photometry, log)[source]

Unpack the data for a single object and correct for Galactic extinction. Also flag pixels which may be affected by emission lines.

fastspecfit.io.unpack_one_stacked_spectrum(iobj, specdata, meta, fphoto, synthphot, ignore_photometry, log)[source]

Unpack the data for a single stacked spectrum. Also flag pixels which may be affected by emission lines.

fastspecfit.io.write_fastspecfit(out, meta, modelspectra=None, outfile=None, specprod=None, coadd_type=None, fphotofile=None, templates=None, emlinesfile=None, fastphot=False, inputz=False, no_smooth_continuum=False, ignore_photometry=False, broadlinefit=True, use_quasarnet=True, constrain_age=False, verbose=True)[source]

Write out.

fastspecfit.mpi

MPI tools.

fastspecfit.mpi._domerge(outfiles, extname='FASTSPEC', survey=None, program=None, outprefix=None, specprod=None, coadd_type=None, mergefile=None, fastphot=False, mp=1)[source]

Support routine to merge a set of input files.

fastspecfit.mpi.merge_fastspecfit(specprod=None, coadd_type=None, survey=None, program=None, healpix=None, tile=None, night=None, sample=None, outsuffix=None, fastphot=False, specprod_dir=None, outdir_data='.', fastfiles_to_merge=None, merge_suffix=None, mergedir=None, supermerge=False, overwrite=False, mp=1)[source]

Merge all the individual catalogs into a single large catalog. Runs only on rank 0.

supermerge - merge previously merged catalogs

fastspecfit.qa

fastspecfit.qa._desiqa_one(args)[source]

Multiprocessing wrapper.

fastspecfit.qa.desiqa_one(data, fastfit, metadata, templates, coadd_type, fphoto, minspecwave=3500.0, maxspecwave=9900.0, minphotwave=0.1, maxphotwave=35.0, emline_snrmin=0.0, nsmoothspec=1, fastphot=False, ignore_photometry=False, stackfit=False, inputz=False, no_smooth_continuum=False, outdir=None, outprefix=None, log=None)[source]

Multiprocessing wrapper to generate QA for a single object.

fastspecfit.qa.fastqa(args=None, comm=None)[source]

Main module.

fastspecfit.qa.parse(options=None)[source]

Parse input arguments.

fastspecfit.qa.qa_fastspec(data, templatecache, fastspec, metadata, coadd_type='healpix', spec_wavelims=(3550, 9900), phot_wavelims=(0.1, 35), fastphot=False, fphoto=None, stackfit=False, outprefix=None, no_smooth_continuum=False, ignore_photometry=False, emline_snrmin=0.0, nsmoothspec=1, outdir=None, log=None, cosmo=None)[source]

QA plot the emission-line spectrum and best-fitting model.

fastspecfit.templates.qa

QA for templates

fastspecfit.templates.qa.qa_bpt(targetclass, fastspecfile=None, png=None)[source]

QA of the fastspec emission-line spectra.

fastspecfit.templates.qa.qa_fastspec_emlinespec(targetclass, fastwave=None, fastflux=None, fastivar=None, fastmeta=None, fastspec=None, fastspecfile=None, CFit=None, EMFit=None, ncol=3, nrow=5, pdffile=None)[source]

QA of the fastspec emission-line spectra.

fastspecfit.templates.qa.qa_fastspec_fullspec(targetclass, fastwave=None, fastflux=None, fastivar=None, fastmeta=None, fastspec=None, fastspecfile=None, CFit=None, EMFit=None, ncol=3, nrow=5, photometric_models=False, pdffile=None)[source]

Full-spectrum QA.

photometric_models - use the fits to the broadband continuum

fastspecfit.templates.qa.qa_parent_sample(samplefile, tilefile, targetclass='lrg', specprod='denali', png=None)[source]

Build QA showing how the parent sample was selected.

fastspecfit.templates.qa.qa_photometry(targetclass, samplefile=None, png_obs=None, png_rest=None, png_rest_bins=None)[source]

QA of the observed- and rest-frame photometry.

fastspecfit.templates.qa.qa_photometry_templates(targetclass, samplefile=None, templatefile=None, ntspace=5, png=None)[source]

Compare the color-color tracks of the templates to the data.

fastspecfit.templates.sample

Code for defining samples and reading the data needed to build templates for the various target classes. Called by bin/desi-templates.

fastspecfit.templates.sample.read_fastspecfit(tilestable, fastspecfit_dir=None, specprod='denali', targetclass='lrg')[source]

Read the fastspecfit output for this production.

fastspecfit.templates.sample.read_parent_sample(samplefile)[source]

Read the output of select_parent_sample

fastspecfit.templates.sample.read_tilestable(tilefile)[source]

Read the output of select_tiles

fastspecfit.templates.sample.select_parent_sample(phot, spec, meta, targetclass='lrg', specprod='denali', deltachi2_cut=40, fastphot_chi2cut=100.0, fastspec_chi2cut=3.0, smoothcorr_cut=10, zobj_minmax=None, absmag_minmax=None, color_minmax=None, return_indices=False, verbose=False, png=None, samplefile=None)[source]

High-level sample selection.

fastspecfit.templates.sample.select_tiles(targetclass, remove_vi=False, min_efftime=None, specprod='denali', outfile=None, png=None)[source]

Select tiles to use. Remove low exposure-time tiles and also remove tiles which are being visually inspected.

/global/cfs/cdirs/desi/sv/vi/TruthTables/Blanc/

fastspecfit.templates.templates

Tools for generating stacked spectra and building templates.

fastspecfit.templates.templates._fastspec_onestack(args)[source]

Multiprocessing wrapper.

fastspecfit.templates.templates._rebuild_fastspec_spectrum(args)[source]

Multiprocessing wrapper.

fastspecfit.templates.templates._spectra_allbins_onetile(args)[source]

Multiprocessing wrapper.

fastspecfit.templates.templates._spectra_onebin_onetile(args)[source]

Multiprocessing wrapper.

fastspecfit.templates.templates._stack_onebin(args)[source]

Multiprocessing wrapper.

fastspecfit.templates.templates.build_templates(fastspecfile, mp=1, minwave=None, maxwave=None, templatefile=None)[source]

Build the final templates.

Called by bin/desi-templates.

fastspecfile is the output of fastspecfit_stacks

fastspecfit.templates.templates.fastspec_onestack(fastmeta, fastspec, flux, ivar, meta, wave, CFit, EMFit, qadir, qaprefix)[source]

Fit one stacked spectrum.

fastspecfit.templates.templates.fastspec_to_desidata(fastspec, wave, flux, ivar)[source]

Convenience wrapper to turn spectra in a fastspecfit-compatible dictionary of DESI data.

fastspecfit.templates.templates.fastspecfit_stacks(stackfile, mp=1, qadir=None, qaprefix=None, fastspecfile=None)[source]

Model stacked spectra using fastspecfit.

Called by bin/desi-templates.

stackfile is the output of stack_in_bins.

fastspecfit.templates.templates.iterative_stack(wave, flux2d, ivar2d=None, constant_ivar=False, maxdiff=0.01, maxiter=500, normwave=4500, smooth=None, debug=False, verbose=False)[source]

Iterative stacking algorithm taken from Bovy, Hogg, & Moustakas 2008. https://arxiv.org/pdf/0805.1200.pdf

fastspecfit.templates.templates.quick_stack(wave, flux2d, ivar2d=None, constant_ivar=False)[source]

Simple inverse-variance-weighted stack.

fastspecfit.templates.templates.rebuild_fastspec_spectrum(fastspec, wave, flux, ivar, CFit, EMFit, full_resolution=False, emlines_only=False, normalize_wave=5500.0, normalize_filter=None, normalize_mag=20.0)[source]

Rebuilding a single fastspec model continuum and emission-line spectrum given a fastspec astropy.table.Table.

full_resolution - returns a full-resolution spectrum in the rest-frame!

fastspecfit.templates.templates.remove_undetected_lines(fastspec, linetable=None, snrmin=3.0, oiidoublet=0.73, siidoublet=1.3, niidoublet=2.936, oiiidoublet=2.875, fastmeta=None, devshift=False)[source]

replace weak or undetected emission lines with their upper limits

devshift - remove individual velocity shifts oiidoublet - [OII] 3726/3729 siidoublet - [SII] 6716/6731 niidoublet - [NII] 6584/6548 oiiidoublet - [OIII] 5007/4959

fastspecfit.templates.templates.spectra_allbins_onetile(fastspecdir, targetclass, tile, bins, minperbin=3)[source]

Find the indices of the objects we care about in a given tile across all the bins of properties.

fastspecfit.templates.templates.spectra_in_bins(tilestable, minperbin=3, targetclass='lrg', specprod='denali', fastspecfit_dir=None, minwave=None, maxwave=None, mp=1, normwave=None, fastphot_in_bins=True, verbose=False, stackfile=None)[source]

Select objects in bins of rest-frame properties.

fastphot_in_bins - also stack the fastphot continuum-fitting results

fastspecfit.templates.templates.spectra_onebin_onetile(fastspecdir, targetclass, tile, bins, minperbin=3, CFit=None, continuumwave=None, index=None, return_index=False)[source]

For Find (and, optionally, read) all the spectra in a given tile and bin of properties.

fastspecfit.templates.templates.stack_in_bins(sample, data, templatewave, mp=1, normwave=4500.0, continuumwave=None, stackfile=None)[source]

Stack spectra in bins given the output of a spectra_in_bins run.

See bin/desi-templates for how to generate the required input.

This routine is basically obsolete – we now use spectra_in_bins directly because there were memory issues passing around the grids of spectra before making the stacks.

fastspecfit.templates.templates.stack_onebin(flux2d, ivar2d, templatewave, normwave, cflux2d, continuumwave)[source]

Build one stack.

fastspecfit.templates.templates.write_binned_stacks(outfile, wave, flux, ivar, metadata=None, cwave=None, cflux=None, fastspec=None, fastmeta=None)[source]

Write out the stacked spectra.

fastspecfit.templates.templates.write_templates(outfile, wave, flux, metadata, weights=None, empca=False)[source]

Write out the final templates

fastspecfit.util

General utilities.

class fastspecfit.util.TabulatedDESI[source]

Class to load tabulated z->E(z) and z->comoving_radial_distance(z) relations within DESI fiducial cosmology (in LSS/data/desi_fiducial_cosmology.dat) and perform the (linear) interpolations at any z.

>>> cosmo = TabulatedDESI()
>>> distance = cosmo.comoving_radial_distance([0.1, 0.2])
>>> efunc = cosmo.efunc(0.3)

The cosmology is defined in https://github.com/abacusorg/AbacusSummit/blob/master/Cosmologies/abacus_cosm000/CLASS.ini and the tabulated file was obtained using https://github.com/adematti/cosmoprimo/blob/main/cosmoprimo/fiducial.py.

Note

Redshift interpolation range is [0, 100].

comoving_radial_distance(z)[source]

Return comoving radial distance, in \(\mathrm{Mpc}/h\).

distance_modulus(z)[source]

Return the distance modulus at the given redshift (Hogg Eq. 24).

efunc(z)[source]

Return \(E(z)\), where the Hubble parameter is defined as \(H(z) = H_{0} E(z)\), unitless.

luminosity_distance(z)[source]

Return luminosity distance, in \(\mathrm{Mpc}/h\).

universe_age(z)[source]

Return the age of the universe at the given redshift.

class fastspecfit.util.ZWarningMask[source]

Mask bit definitions for zwarning.

WARNING on the warnings: not all of these are implemented yet.

#- TODO: Consider using something like desispec.maskbits to provide a more #- convenient wrapper class (probably copy it here; don’t make a dependency) #- That class as-is would bring in a yaml dependency.

fastspecfit.util._trapz_rebin(x, y, edges, results)

Numba-friendly version of trapezoidal rebinning

See redrock.rebin.trapz_rebin() for input descriptions. results is pre-allocated array of length len(edges)-1 to keep results

fastspecfit.util.air2vac(airwave)[source]

http://www.astro.uu.se/valdwiki/Air-to-vacuum%20conversion

fastspecfit.util.centers2edges(centers)[source]

Convert bin centers to bin edges, guessing at what you probably meant

Parameters:

centers (array) – bin centers,

Returns:

bin edges, lenth = len(centers) + 1

Return type:

array

fastspecfit.util.find_minima(x)[source]

Return indices of local minima of x, including edges.

The indices are sorted small to large.

Note

this is somewhat conservative in the case of repeated values: find_minima([1,1,1,2,2,2]) -> [0,1,2,4,5]

Parameters:

x (array-like) – The data array.

Returns:

The indices.

Return type:

(array)

fastspecfit.util.ivar2var(ivar, clip=0.001, sigma=False, allmasked_ok=False)[source]

Safely convert an inverse variance to a variance. Note that we clip at 1e-3 by default, not zero.

fastspecfit.util.minfit(x, y, return_coeff=False)[source]

Fits y = y0 + ((x-x0)/xerr)**2

See redrock.zwarning.ZWarningMask.BAD_MINFIT for zwarn failure flags

Parameters:
  • x (array) – x values.

  • y (array) – y values.

Returns:

(x0, xerr, y0, zwarn) where zwarn=0 is good fit.

Return type:

(tuple)

fastspecfit.util.mwdust_transmission(ebv, filtername)[source]

Convert SFD E(B-V) value to dust transmission 0-1 given the bandpass.

Parameters:
  • ebv (float or array-like) – SFD E(B-V) value(s)

  • filtername (str) – Filter name, e.g., ‘decam2014-r’.

Returns:

Scalar or array (same as ebv input), Milky Way dust transmission 0-1.

Note

This function tabulates the total-to-selective extinction ratio, k_X=A(X)/E(B-V) for many different filter bandpasses, X, where A(X)=-2.5*log10(transmission in X band). And so given the ebv, it returns mwdust_transmission=10**(-0.4*k_X*ebv).

Returns:

scalar, total extinction A(band) = -2.5*log10(transmission(band))

Notes

Based on desiutil.dust.mwdust_transmission.

fastspecfit.util.trapz_rebin(x, y, xnew=None, edges=None)[source]

Rebin y(x) flux density using trapezoidal integration between bin edges

Notes

y is interpreted as a density, as is the output, e.g.

>>> x = np.arange(10)
>>> y = np.ones(10)
>>> trapz_rebin(x, y, edges=[0,2,4,6,8])  #- density still 1, not 2
array([ 1.,  1.,  1.,  1.])
Parameters:
  • x (array) – input x values.

  • y (array) – input y values.

  • edges (array) – (optional) new bin edges.

Returns:

integrated results with len(results) = len(edges)-1

Return type:

array

Raises:

ValueError – if edges are outside the range of x or if len(x) != len(y)