FastSpecFit API¶
fastspecfit¶
Tools for fast stellar continuum, emission-line, and broadband photometric modeling.
fastspecfit.continuum¶
Methods and tools for continuum-fitting.
- class fastspecfit.continuum.ContinuumFit(ssptemplates=None, metallicity='Z0.0190', minwave=None, maxwave=300000.0, nolegend=False, cache_vdisp=True, solve_vdisp=False, cache_SSPgrid=True, constrain_age=True, mapdir=None)[source]¶
- _fnnls_parallel(modelflux, flux, ivar, xparam=None, debug=False, interpolate_coeff=False, xlabel=None)[source]¶
Wrapper on fnnls to set up the multiprocessing.
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_fastphot(data, minbands=3)[source]¶
Fit the broad photometry.
- Parameters
- Returns
astropy.table.Table
– Table with all the continuum-fitting results with columns documented in init_phot_output... note:: – See https://github.com/mikeiovine/fast-nnls for the fNNLS algorithm(s).
- continuum_specfit(data)[source]¶
Fit the non-negative stellar continuum of a single spectrum.
- Parameters
data (
dict
) – Dictionary of input spectroscopy (plus ancillary data) populated byfastspecfit.io.DESISpectra.read_and_unpack()
.- Returns
Table with all the continuum-fitting results with columns documented in
self.init_spec_output()
.- Return type
Notes
Consider using cross-correlation to update the redrock redshift.
We solve for velocity dispersion if solve_vdisp=True or ((SNR_B>3 or SNR_R>3) and REDSHIFT<1).
- kcorr_and_absmag(data, continuum, coeff, snrmin=2.0)[source]¶
Computer K-corrections, absolute magnitudes, and a simple stellar mass.
- class fastspecfit.continuum.ContinuumTools(ssptemplates=None, metallicity='Z0.0190', minwave=None, maxwave=300000.0, sspversion='v1.0', mapdir=None)[source]¶
Tools for dealing with stellar continua.
- Parameters
metallicity (
str
, optional, defaults to Z0.0190.) – Stellar metallicity of the SSPs. Currently fixed at solar metallicity, Z=0.0190.minwave (
float
, optional, defaults to None) – Minimum SSP wavelength to read into memory. IfNone
, the minimum available wavelength is used (around 100 Angstrom).maxwave (
float
, optional, defaults to 6e4) – Maximum SSP wavelength to read into memory.note:: (..) – Need to document all the attributes.
- SSP2data(_sspflux, _sspwave, redshift=0.0, AV=None, vdisp=None, cameras=['b', 'r', 'z'], specwave=None, specres=None, coeff=None, south=True, synthphot=True, debug=False)[source]¶
Workhorse routine to turn input SSPs 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 dust reddening - apply velocity dispersion broadening - apply the resolution matrix - synthesize photometry
It also naturally handles SSPs which have been precomputed on a grid of reddening or 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?).
- build_linemask(wave, flux, ivar, redshift=0.0, nsig=7.0)[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.
- convolve_vdisp(sspflux, vdisp)[source]¶
Convolve by the velocity dispersion.
- Parameters
sspflux –
vdisp –
Notes
- dust_attenuation(wave, AV, test=False)[source]¶
Compute the dust attenuation curve A(lambda)/A(V) from Charlot & Fall 2000.
- ToDo: add a UV bump and IGM attenuation!
https://gitlab.lam.fr/cigale/cigale/-/blob/master/pcigale/sed_modules/dustatt_powerlaw.py#L42
- estimate_linesigmas(wave, flux, ivar, redshift=0.0, png=None, refit=True)[source]¶
Estimate the velocity width from potentially strong, isolated lines.
- static get_dn4000(wave, flam, flam_ivar=None, redshift=None, rest=True)[source]¶
Compute DN(4000) and, optionally, the inverse variance.
- Parameters
wave –
flam –
flam_ivar –
redshift –
rest –
Notes
If rest`=``False` then redshift input is required.
Require full wavelength coverage over the definition of the index.
See eq. 11 in Bruzual 1983 (https://articles.adsabs.harvard.edu/pdf/1983ApJ…273..105B) but with the “narrow” definition of Balogh et al. 1999.
- static parse_photometry(bands, maggies, lambda_eff, ivarmaggies=None, nanomaggies=True, nsigma=2.0, min_uncertainty=None, debug=False)[source]¶
Parse input (nano)maggies to various outputs and pack into a table.
- Parameters
erg/s/cm2/A (flam - 10-17) –
erg/s/cm2/Hz (fnu - 10-17) –
mag (abmag - AB) –
maggies (nanomaggies - input maggies are actually 1e-9) –
limit (nsigma - magnitude) –
- Return type
phot - photometric table
Notes
- smooth_and_resample(sspflux, sspwave, specwave=None, specres=None)[source]¶
Given a single template, apply the resolution matrix and resample in wavelength.
- Parameters
sspflux (
numpy.ndarray
[npix]) – Input (model) spectrum.sspwave (
numpy.ndarray
[npix]) – Wavelength array corresponding to sspflux.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.vdisp (
float
, optional, defaults to None) – Velocity dispersion broadening factor [km/s].pixkms (
float
, optional, defaults to None) – Pixel size of input spectra [km/s].
- 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.
- smooth_continuum(wave, flux, ivar, redshift, medbin=150, smooth_window=50, smooth_step=10, maskkms_uv=3000.0, maskkms_balmer=1000.0, maskkms_narrow=200.0, linemask=None, png=None)[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 typebool
, 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._fnnls_continuum(myargs)[source]¶
Multiprocessing wrapper for fnnls_continuum.
- fastspecfit.continuum.fnnls_continuum(ZZ, xx, flux=None, ivar=None, modelflux=None, get_chi2=False)[source]¶
Fit a stellar continuum using fNNLS.
- Parameters
- Returns
Notes
This function is a simple wrapper on fastspecfit.fnnls.fnnls(); see the ContinuumFit.fnnls_continuum method for documentation.
The arguments
flux
,ivar
andmodelflux
are only used whenget_chi2=True
.
fastspecfit.emlines¶
Methods and tools for fitting emission lines.
- class fastspecfit.emlines.EMLineFit(nball=10, chi2_default=0.0, minwave=3000.0, maxwave=10000.0, pixkms=10.0, nolegend=False, ssptemplates=None, mapdir=None)[source]¶
Class to fit an emission-line spectrum.
https://docs.astropy.org/en/stable/modeling/example-fitting-constraints.html#tied
https://docs.astropy.org/en/stable/modeling/compound-models.html#parameters
- chi2(bestfit, emlinewave, emlineflux, emlineivar, continuum_model=None, return_dof=False)[source]¶
Compute the reduced chi^2.
- emlinemodel_bestfit(specwave, specres, fastspecfit_table, redshift=None)[source]¶
Wrapper function to get the best-fitting emission-line model from an fastspecfit table (to be used to build QA).
- fit(data, continuummodel, smooth_continuum, synthphot=True, maxiter=5000, accuracy=0.01, verbose=False)[source]¶
Perform the fit minimization / chi2 minimization.
EMLineModel object FC - ContinuumFit object
ToDo: need to take into account the instrumental velocity width when computing integrated fluxes…
- class fastspecfit.emlines.EMLineModel(*args: Any, **kwargs: Any)[source]¶
Class to model the emission-line spectra.
- _emline_spectrum(*lineargs)[source]¶
Simple wrapper to build an emission-line spectrum.
# build the emission-line model [erg/s/cm2/A, observed frame]
- class fastspecfit.emlines.FastLevMarLSQFitter(*args: Any, **kwargs: Any)[source]¶
- Adopted in part from
https://github.com/astropy/astropy/blob/v4.2.x/astropy/modeling/fitting.py#L1580-L1671
Copyright (c) 2011-2021, Astropy Developers
See https://github.com/astropy/astropy/issues/12089 for details.
- fitter_to_model_params(model, fps)[source]¶
Constructs the full list of model parameters from the fitted and constrained parameters.
- model_to_fit_params(model)[source]¶
Convert a model instance’s parameter array to an array that can be used with a fitter that doesn’t natively support fixed or tied parameters. In particular, it removes fixed/tied parameters from the parameter array. These may be a subset of the model parameters, if some of them are held constant or tied.
- fastspecfit.emlines._tie_lines(model)[source]¶
Tie emission lines together using sensible constrains based on galaxy & QSO physics.
- fastspecfit.emlines._tie_nii_amp(model)[source]¶
[N2] (4–>2): airwave: 6548.0488 vacwave: 6549.8578 emissivity: 2.022e-21 [N2] (4–>3): airwave: 6583.4511 vacwave: 6585.2696 emissivity: 5.949e-21
- fastspecfit.emlines._tie_oii_red_amp(model)[source]¶
[O2] (4–>2): airwave: 7319.9849 vacwave: 7322.0018 emissivity: 3.229e-22 [O2] (4–>3): airwave: 7330.7308 vacwave: 7332.7506 emissivity: 1.692e-22
fastspecfit.fastspecfit¶
See sandbox/running-fastspecfit for examples.
- fastspecfit.fastspecfit._assign_units_to_columns(fastfit, metadata, Spec, CFit, EMFit=None, fastphot=False)[source]¶
Assign astropy units to output tables.
- fastspecfit.fastspecfit.desiqa_one(CFit, EMFit, data, fastfit, metadata, coadd_type, fastphot=False, outdir=None, outprefix=None)[source]¶
Multiprocessing wrapper to generate QA for a single object.
- 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
orNone
) – 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.fastphot_one(iobj, data, out, meta, CFit)[source]¶
Multiprocessing wrapper to run
fastphot()
on a single object.
- fastspecfit.fastspecfit.fastspec(args=None, comm=None)[source]¶
Main fastspec script.
This script is the engine to model one or more DESI spectra. It initializes the
ContinuumFit
andEMLineFit
classes, 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
orNone
) – 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, CFit, EMFit, verbose=False)[source]¶
Multiprocessing wrapper to run
fastspec()
on a single object.
fastspecfit.fnnls¶
https://github.com/mikeiovine/fast-nnls.git
This is a Python implementation of the algorithm described in the paper “A Fast Non-Negativity-Constrained Least Squares Algorithm” by Rasmus Bro and Sumen De Jong.
Give a matrix A and a vector y, this algorithm solves argmin_x ||Ax - y||.
At the time of writing this, there are no readily available Python bindings for this algorithm that I know of. scipy.optimize.nnls implements the slower Lawson-Hanson version.
## Usage
A = np.array([[1, 2], [3, 4], [5, 6]]) y = np.array([1, 2, 3]) AtA = A.T.dot(A) Aty = A.T.dot(y) x = fnnls(AtA, Aty)
- fastspecfit.fnnls.fnnls(AtA, Aty, epsilon=None, iter_max=None)[source]¶
Given a matrix A and vector y, find x which minimizes the objective function f(x) = ||Ax - y||^2. This algorithm is similar to the widespread Lawson-Hanson method, but implements the optimizations described in the paper “A Fast Non-Negativity-Constrained Least Squares Algorithm” by Rasmus Bro and Sumen De Jong.
Note that the inputs are not A and y, but are A^T * A and A^T * y
This is to avoid incurring the overhead of computing these products many times in cases where we need to call this routine many times.
- Parameters
AtA (numpy.ndarray) – A^T * A. See above for definitions. If A is an (m x n) matrix, this should be an (n x n) matrix.
Aty (numpy.ndarray) – A^T * y. See above for definitions. If A is an (m x n) matrix and y is an m dimensional vector, this should be an n dimensional vector.
epsilon (float) – Anything less than this value is consider 0 in the code. Use this to prevent issues with floating point precision. Defaults to the machine precision for doubles.
iter_max (int, optional) – Maximum number of inner loop iterations. Defaults to 30 * [number of cols in A] (the same value that is used in the publication this algorithm comes from).
fastspecfit.html¶
Code to support building HTML visualizations.
fastspecfit.io¶
Tools for reading DESI spectra and reading and writing fastspecfit files.
- fastspecfit.io.read_fastspecfit(fastfitfile, rows=None, columns=None)[source]¶
Read the fitting results.
fastspecfit.mpi¶
MPI tools.
- fastspecfit.mpi.backup_logs(logfile)[source]¶
Move logfile -> logfile.0 or logfile.1 or logfile.n as needed
TODO: make robust against logfile.abc also existing
- fastspecfit.mpi.group_redrockfiles(specfiles, maxnodes=256, comm=None, makeqa=False)[source]¶
Group redrockfiles to balance runtimes
- Parameters
specfiles – list of spectra filepaths
- Options:
maxnodes: split the spectra into this number of nodes comm: MPI communicator
- Returns (groups, ntargets, grouptimes):
groups: list of lists of indices to specfiles
list of number of targets per group
grouptimes: list of expected runtimes for each group
- fastspecfit.mpi.merge_fastspecfit(specprod=None, coadd_type=None, survey=None, program=None, healpix=None, tile=None, night=None, outsuffix=None, fastphot=False, specprod_dir=None, outdir_data='.', 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.mpi.weighted_partition(weights, n)[source]¶
Partition
weights
inton
groups with approximately same sum(weights).- Parameters
weights – array-like weights
n – number of groups
- Returns (groups, groupweights):
groups: list of lists of indices of weights for each group groupweights: sum of weights assigned to each group
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.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.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.templates¶
Tools for generating stacked spectra and building templates.
- 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.util¶
General utilities.
- 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.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)[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.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)