API Documentation¶
lightcurve_fitting.bolometric module¶
- lightcurve_fitting.bolometric.blackbody_lstsq(epoch1, z, p0=None, T_range=(1.0, 100.0), R_range=(0.01, 1000.0), cutoff_freq=inf)[source]¶
Fit a blackbody spectrum to a spectral energy distribution using \(χ^2\) minimization
- Parameters:
- epoch1lightcurve_fitting.lightcurve.LC
A single “epoch” of photometry that defines the observed spectral energy distribution
- zfloat
Redshift between the blackbody (rest frame) and the filters (observed frame)
- p0list, tuple, array-like, optional
Initial guess for [temperature (kK), radius (1000 Rsun)]. Default:
[10., 10.]
- T_rangetuple, list, array-like, optional
Range of allowed temperatures (in kilokelvins) in the prior. Default:
(1., 100.)
- R_rangetuple, list, array-like, optional
Range of allowed radii (in 1000 solar radii) in the prior. Default:
(0.01, 1000.)
- cutoff_freqfloat, optional
Cutoff frequency (in terahertz) for a modified blackbody spectrum (see https://doi.org/10.3847/1538-4357/aa9334)
- Returns:
- tempfloat
Best-fit blackbody temperature in kilokelvins
- radiusfloat
Best-fit blackbody radius in units of 1000 solar radii
- dtempfloat
Uncertainty in the temperature in kilokelvins
- dradfloat
Uncertainty in the radius in units of 1000 solar radii
- lumfloat
Blackbody luminosity implied by the best-fit temperature and radius
- dlumfloat
Uncertainty in the blackbody luminosity
- L_optfloat
Pseudobolometric luminosity from integrating the best-fit blackbody spectrum over the U-I filters
- lightcurve_fitting.bolometric.calc_colors(epoch1, colors)[source]¶
Calculate colors from a spectral energy distribution
- Parameters:
- epoch1lightcurve_fitting.lightcurve.LC
A single “epoch” of photometry that defines the observed spectral energy distribution
- colorslist
A list of colors to calculate, e.g.,
['U-B', 'B-V', 'g-r', 'r-i']
- Returns:
- magslist
A list of the calculated colors
- dmagslist
A list of uncertainties on the calculated colors, from propagating errors on the photometry points
- lolimslist
A list of booleans: True if the first filter is a nondetection or is absent from
epoch1
- uplimslist
A list of booleans: True if the second filter is a nondetection or is absent from
epoch1
- lightcurve_fitting.bolometric.calculate_bolometric(lc, z=0.0, outpath='.', res=1.0, nwalkers=10, burnin_steps=200, steps=100, priors=None, save_table_as=None, min_nfilt=3, cutoff_freq=inf, show=False, colors=None, do_mcmc=True, save_chains=False, use_sigma=False, sigma_type='relative', also_group_by=())[source]¶
Calculate the full bolometric light curve from a table of broadband photometry.
To manually specify epochs for some or all points, add a unique number in the ‘epoch’ column of the
lc
.- Parameters:
- lclightcurve_fitting.lightcurve.LC
Table of broadband photometry including columns “MJD”, “mag”, “dmag”, “filter”
- zfloat, optional
DEPRECATED: Include the redshift in lc.meta[‘redshift’] instead.
- outpathstr, optional
Directory to which to save the corner plots and MCMC chains. Default: current directory
- resfloat, optional
Approximate resolution for grouping, in days. Default: 1 day.
- nwalkersint, optional
Number of walkers (chains) to use for fitting. Default: 10
- burnin_stepsint, optional
Number of MCMC steps before convergence. This part of the history is discarded. Default: 200
- stepsint, optional
Number of MCMC steps after convergence. This part of the history is used to calculate paramers. Default: 100.
- priorslist, optional
Prior probability distributions for temperature (in kilokelvins) and radius (in 1000 solar radii). Available priors:
models.UniformPrior
,models.LogUniformPrior
,models.GaussianPrior
. Default:T = Uniform(1., 100.)
,R = LogUniform(0.01, 1000.)
,σ = Gaussian(0., 10., stddev=1.)
.- save_table_asstr, optional
Filename to which to save the output table of blackbody parameters and bolometric luminosities
- min_nfiltint, optional
Minimum number of distinct observed filters required to calculate a luminosity. Default: 3
- cutoff_freqfloat, optional
Cutoff frequency (in terahertz) for a modified blackbody spectrum (see https://doi.org/10.3847/1538-4357/aa9334)
- showbool, optional
Plot the chain histories, and display all plots at the end of the script. Default: only save the corner plot
- colorslist, optional
A list of colors to calculate, e.g.,
['U-B', 'B-V', 'g-r', 'r-i']
- do_mcmcbool, optional
If True (default), also fit the spectral energy distribution with an MCMC routine. This is slower but gives more realistic uncertainties.
- save_chainsbool, optional
If True, save the MCMC chain histories to the directory
outpath
. Default: only save the corner plot- use_sigmabool, optional
Include an intrinsic scatter parameter in the MCMC fit. Default: False.
- sigma_typestr, optional
If ‘relative’ (default), sigma will be in units of the individual photometric uncertainties. If ‘absolute’, sigma will be in units of the median photometric uncertainty.
- also_group_bylist, optional
Group by these columns in addition to epoch
- Returns:
- t0lightcurve_fitting.lightcurve.LC
Table containing the blackbody parameters, bolometric luminosities, and (optionally) colors
- lightcurve_fitting.bolometric.group_by_epoch(lc, res=1.0, also_group_by=())[source]¶
Group a light curve into epochs that will be treated as single spectral energy distributions.
To manually specify epochs for some or all points, add a unique number in the ‘epoch’ column of the
lc
.- Parameters:
- lclightcurve_fitting.lightcurve.LC
Table containing the observed photometry. Must contain times in an “MJD” column.
- resfloat, optional
Approximate resolution for grouping, in days. Default: 1 day.
- also_group_bylist, tuple, optional
Column names to group by, in addition to epoch. Default: ().
- Returns:
- groupslist of lightcurve_fitting.lightcurve.LC
Iterable containing single-epoch spectral energy distributions, sorted by average MJD
- lightcurve_fitting.bolometric.integrate_sed(epoch1)[source]¶
Directly integrate a spectral energy distribution using the trapezoidal rule
- Parameters:
- epoch1lightcurve_fitting.lightcurve.LC
A single “epoch” of photometry that defines the observed spectral energy distribution
- Returns:
- L_intfloat
Luminosity in watts
- lightcurve_fitting.bolometric.median_and_unc(x, perc_contained=68.0)[source]¶
Calculate the equal-tailed credible interval, centered on the median, for a sample of data
- Parameters:
- xarray-like
The data sample
- perc_containedfloat, optional
The percentage of the probability to be contained in the interval. Default: 68% (1σ)
- Returns:
- medianfloat
The median of
x
- lowerfloat
The lower boundary of the credible interval
- upperfloat
The upper boundary of the credible interval
- lightcurve_fitting.bolometric.plot_bolometric_results(t0, save_plot_as=None, xcol=None, log=False)[source]¶
Plot the parameters and bolometric light curves that result from fitting the spectral energy distribution
- Parameters:
- t0lightcurve_fitting.lightcurve.LC
Table containing the results from fitting the spectral energy distribution
- save_plot_asstr, optional
Filename to which to save the plot.
- xcolstr, optional
Column to plot on the x-axis. Default:
'phase'
ift0.meta['refmjd']
andt0.meta['redshift']
are given, otherwise'MJD'
. (See the docs forLC.calcPhase()
.)- logbool, optional
Use a logarithmic x-axis. Only useful if
xcol='phase'
. Default: False.
- Returns:
- figmatplotlib.pyplot.Figure
Figure object containing the plot
- lightcurve_fitting.bolometric.plot_chain(chain, labels=None)[source]¶
Plot the chain histories for an MCMC run
- Parameters:
- chainarray-like
Array containing the MCMC chain history, e.g., from
emcee.EnsembleSampler.chain()
- labelsiterable, optional
A list of axis labels for each parameter in the chain
- Returns:
- figmatplotlib.pyplot.Figure
Figure object containing the chain history plots
- lightcurve_fitting.bolometric.plot_color_curves(t, colors=None, fmt='o', limit_length=0.1, xcol='MJD')[source]¶
Plot the color curves calculated by
calculate_bolometric()
.- Parameters:
- tlightcurve_fitting.lightcurve.LC
Output table from
calculate_bolometric()
- colorslist, optional
List of colors to plot, e.g.,
['g-r', 'r-i']
. By default, plot all recognizable colors.- fmtstr, optional
Format string, passed to
matplotlib.pyplot.errorbar()
- limit_lengthfloat, optional
Length (in data units) of the upper and lower limit markers
- Returns:
- figmatplotlib.pyplot.Figure
Figure object containing the plot
- lightcurve_fitting.bolometric.pseudo(temp, radius, z, filter0=<filter I>, filter1=<filter U>, cutoff_freq=inf)[source]¶
Integrate a blackbody spectrum between two broadband filters to produce a pseudobolometric luminosity
- Parameters:
- tempfloat, array-like
Blackbody temperatures in kilokelvins
- radiusfloat, array-like
Blackbody radii in units of 1000 solar radii
- zfloat
Redshift between the blackbody (rest frame) and the filters (observed frame)
- filter0, filter1lightcurve_fitting.filters.Filter, optional
Filters to integrate between, where
filter1
is bluer thanfilter0
. Default: U to I.- cutoff_freqfloat, optional
Cutoff frequency (in terahertz) for a modified blackbody spectrum (see https://doi.org/10.3847/1538-4357/aa9334)
- Returns:
- L_optfloat, array-like
Pseudobolometric luminosity in watts
- lightcurve_fitting.bolometric.spectrum_corner(spectrum, epoch1, sampler_flatchain, z=0.0, ebv=0.0, spectrum_kwargs=None, use_sigma=False, labels=None, freq_min=100.0, freq_max=1000.0, save_plot_as='')[source]¶
Plot the posterior distributions in a corner (pair) plot, with an inset showing the observed and model SEDs.
- Parameters:
- spectrumfunction
Function describing the spectrum. The first argument must be frequency in THz, and it must return spectral luminosity in watts per hertz. For a blackbody, you can use
models.planck_fast()
.- epoch1lightcurve_fitting.lightcurve.LC
A single “epoch” of photometry that defines the observed spectral energy distribution
- sampler_flatchainarray-like
2D array containing the aggregated MCMC chain histories
- zfloat, optional
Redshift between the emission source and the observed filter. Default: 0.
- ebvfloat, array-like, optional
Selective extinction \(E(B-V)\) in magnitudes, evaluated using a Fitzpatrick (1999) extinction law with \(R_V=3.1\). Its shape must be broadcastable to any array-like arguments. Default: 0.
- spectrum_kwargsdict, optional
Keyword arguments to be passed to the
spectrum
function.- use_sigmabool, optional
Include an intrinsic scatter parameter. Default: False.
- labelslist, optional
Axis labels for the chain histories and corner plot.
- freq_min, freq_maxfloat, optional
Minimum and maximum frequencies for the SED panel in the corner plot. Default: (100., 1000.) or observed range.
- save_plot_asstr, optional
Filename to which to save the resulting plot
- Returns:
- f4matplotlib.pyplot.Figure
Figure containing the corner and SED plots.
- lightcurve_fitting.bolometric.spectrum_mcmc(spectrum, epoch1, priors, starting_guesses, z=0.0, ebv=0.0, spectrum_kwargs=None, show=False, outpath='.', nwalkers=10, burnin_steps=200, steps=100, save_chains=False, use_sigma=False, sigma_type='relative', labels=None, freq_min=100.0, freq_max=1000.0)[source]¶
Fit the given spectral energy distribution to an epoch of photometry using a Markov-chain Monte Carlo routine
- Parameters:
- spectrumfunction
Function describing the spectrum. The first argument must be frequency in THz, and it must return spectral luminosity in watts per hertz. For a blackbody, you can use
models.planck_fast()
.- epoch1lightcurve_fitting.lightcurve.LC
A single “epoch” of photometry that defines the observed spectral energy distribution
- priorslist
Prior probability distributions for each model parameter (and sigma, if used). Available priors:
models.UniformPrior
(default),models.LogUniformPrior
,models.GaussianPrior
- starting_guessesarray-like
Initial guesses for each input parameter to
spectrum
. Must have shape (nwalkers, nparameters).- zfloat, optional
Redshift between the emission source and the observed filter. Default: 0.
- ebvfloat, array-like, optional
Selective extinction \(E(B-V)\) in magnitudes, evaluated using a Fitzpatrick (1999) extinction law with \(R_V=3.1\). Its shape must be broadcastable to any array-like arguments. Default: 0.
- spectrum_kwargsdict, optional
Keyword arguments to be passed to the
spectrum
function.- showbool, optional
Plot the chain histories, and display all plots at the end of the script. Default: only save the corner plot
- outpathstr, optional
Directory to which to save the corner plots. Default: current directory
- nwalkersint, optional
Number of walkers (chains) to use for fitting. Default: 10
- burnin_stepsint, optional
Number of MCMC steps before convergence. This part of the history is discarded. Default: 200
- stepsint, optional
Number of MCMC steps after convergence. This part of the history is used to calculate parameters. Default: 100.
- save_chainsbool, optional
If True, save the MCMC chain histories to the directory
outpath
. Default: only save the corner plot- use_sigmabool, optional
Include an intrinsic scatter parameter. Default: False.
- sigma_typestr, optional
If ‘relative’ (default), sigma will be in units of the individual photometric uncertainties. If ‘absolute’, sigma will be in units of the median photometric uncertainty.
- labelslist, optional
Axis labels for the chain histories and corner plot.
- freq_min, freq_maxfloat, optional
Minimum and maximum frequencies for the SED panel in the corner plot. Default: (100., 1000.) or observed range.
- Returns:
- sampleremcee.EnsembleSampler
Sampler object containing the results of the MCMC fit
- lightcurve_fitting.bolometric.stefan_boltzmann(temp, radius, dtemp=None, drad=None, covTR=None)[source]¶
Calculate blackbody luminosity and associated uncertainty using the Stefan-Boltzmann law
- Parameters:
- tempfloat, array-like
Temperature in kilokelvins.
- radiusfloat, array-like
Radius in units of 1000 solar radii.
- dtempfloat, array-like, optional
Uncertainty in the temperature in kilokelvins
- dradfloat, array-like, optional
Uncertainty in the radius in units of 1000 solar radii
- covTRfloat, array-like, optional
Covariance between the temperature and radius
- Returns:
- lumfloat, array-like
Luminosity in watts
- dlumfloat, array-like
Uncertainty in the luminosity in watts (returned only if dtemp, drad, and covTR are given)
lightcurve_fitting.filters module¶
- class lightcurve_fitting.filters.Filter(names, color='k', offset=0, system=None, fnu=3.631e-23, filename='', angstrom=False, linecolor=None, textcolor=None, mec=None, italics=True)[source]¶
Bases:
object
A broadband photometric filter described by its transmission function and its associated photometric system
- Parameters:
- namesstr, list
One or more names for the filter. The first name is used by default but other names are recognized.
- colorstr, tuple, optional
The color used when plotting photometry in this filter. Default: black.
- offsetfloat, optional
When plotting, offset photometry in this filter by this amount (in magnitudes if plotting magnitudes)
- systemstr, optional
Photometric system. Only used for grouping filters in legends.
- fnufloat, optional
Zero-point flux for magnitudes in this filter, in W/(m^2 Hz). If
fnu = None
andsystem in ['Gunn', 'ATLAS', 'Gaia', 'MOSFiT']
, assume the AB zero point. Otherwise iffnu = None
, converting to flux will raise an error.- filenamestr, optional
Filename containing the transmission function. If none is supplied, converting to flux will raise an error.
- angstrombool, optional
If False (default), the transmission function wavelengths are assumed to be in nanometers. If True, they are given in ångströms.
- linecolorstr, tuple, optional
The line color used when plotting photometry in this filter. Default: same as
color
.- textcolorstr, tuple, optional
The color used when printing the name of this filter. Default: same as
linecolor
.- mecstr, tuple, optional
The marker edge color used when plotting photometry in this filter. Default: same as
linecolor
.- italicsbool, optional
Italicize the filter name when used with LaTeX. Default: True.
- Attributes:
- namestr
The default name of the filter
- nameslist
A list of aliases for the filter
- charstr
A single-character identifier for the filter
- colorstr, tuple
The color used when plotting photometry in this filter
- linecolorstr, tuple
The line and marker edge color used when plotting photometry in this filter
- textcolorstr, tuple
The color used when printing the name of this filter
- italicsbool
Italicize the filter name when used with LaTeX
- mecstr, tuple
The marker edge color used when plotting photometry in this filter
- plotstyledict
Contains all keyword arguments for plotting photometry in this filter
- offsetfloat
When plotting, offset photometry in this filter by this amount (in magnitudes if plotting magnitudes)
- systemstr
The photometric system for magnitudes in this filter
- fnufloat
The zero-point flux for magnitudes in this filter, in watts per square meter per hertz
- m0float
The zero point for magnitudes in this filter, i.e., \(m_0 = 2.5 \log_{10}(F_ν)\)
- M0float
The zero point for absolute magnitudes in this filter, i.e., \(M_0 = m_0 + 90.19\)
- filenamestr
Filename containing the transmission function
- angstrombool
If False (default), the transmission function wavelengths are assumed to be in nanometers. If True, they are given in ångströms.
- transastropy.table.Table
The transmission function
- freq_efffloat
The effective frequency of the filter (in terahertz), i.e., \(ν_0 = \frac{1}{Δν} \int ν T(ν) dν\)
- dfreqfloat
The effective bandwidth of the filter (in terahertz), i.e., \(Δν = \int T(ν) dν\)
- freq_rangetuple
The approximate edges of the filter transmission function (in terahertz), i.e., \(ν_0 \pm Δν\)
- property dfreq¶
- property dwl¶
- extinction(ebv, rv=3.1, z=0.0)[source]¶
Extinction \(A_\lambda\) at the effective wavelength of this filter
- Parameters:
- ebvarray-like
Selective extinction \(E(B-V)\) in magnitudes
- rvfloat, optional
Ratio of total to selective extinction \(R_V\). Default: 3.1.
- zfloat, optional
Redshift between the dust and the observed filter. Default: 0 (appropriate for Milky Way extinction).
- Returns:
- extinctionfloat
Extinction at the effective wavelength of this filter in magnitudes
- property freq_eff¶
- property freq_range¶
- order = ['FUV', 'NUV', 'UVW2', 'UVM2', 'UVW1', 'u', 'U_S', 'U', 'B', 'B_S', 'g', 'g-DECam', 'c', 'V', 'V_S', 'Itagaki', 'white', 'unfilt.', 'G', 'Kepler', 'TESS', 'DLT40', 'w', 'o', 'r', 'r-DECam', 'R', 'i', 'i-DECam', 'I', 'z_s', 'z', 'z-DECam', 'y', 'y-DECam', 'J', 'H', 'K', 'L', 'F070W', 'F090W', 'F115W', 'F150W', 'F182M', 'F200W', 'F250M', 'F277W', 'F300M', 'F335M', 'F356W', 'F360M', 'F444W', 'F560W', 'F770W', 'F1000W', 'F1130W', 'F1280W', 'F1500W', 'F1800W', 'F2100W', 'F2550W', 'pseudobolometric, curve_fit', 'pseudobolometric, MCMC', 'pseudobolometric, integration', 'bolometric, curve_fit', 'bolometric, MCMC', 'unknown']¶
The names of recognized filters listed in (approximate) decreasing order of effective frequency
- read_curve(show=False, force=False)[source]¶
Read the transmission function from
self.filename
and store it inself.trans
- Parameters:
- showbool, optional
If True, also plot the transmission function
- forcebool, optional
If True, reread the transmission function from
self.filename
even if is already stored inself.trans
- spectrum(freq, lum, z=0.0, ebv=0.0)[source]¶
Return the average value of the given spectrum in this filter
The limits of integration come from the input spectrum, so if the spectrum does not cover the full filter, this function returns the average Lnu within the overlapping region, i.e., it does not extrapolate.
- Parameters:
- freqfloat or array-like
Frequency of the input spectrum in THz
- lumfloat or array-like
Spectral luminosity (Lnu) or flux (Fnu) of the input spectrum
- zfloat, optional
Redshift between the spectrum and the filter. Default: 0.
- ebvfloat, array-like, optional
Selective extinction \(E(B-V)\) in magnitudes, evaluated using a Fitzpatrick (1999) extinction law with \(R_V=3.1\). Default: 0.
- Returns:
- Lnufloat or array-like
Average spectral luminosity (Lnu) or flux (Fnu) in the filter
- synthesize(spectrum, *args, z=0.0, ebv=0.0, **kwargs)[source]¶
Returns the average Lnu of the given spectrum in this filter
- Parameters:
- spectrumfunction
Function describing the spectrum. The first argument must be frequency in THz, and it must return spectral luminosity in watts per hertz. All arguments are passed to this function except for keywords z and ebv.
- zfloat, optional
Redshift between the emission source and the observed filter. Default: 0.
- ebvfloat, array-like, optional
Selective extinction \(E(B-V)\) in magnitudes, evaluated using a Fitzpatrick (1999) extinction law with \(R_V=3.1\). Its shape must be broadcastable to any array-like arguments. Default: 0.
- Returns:
- Lnufloat or array-like
Average spectral luminosity in the filter in watts per hertz
- property trans¶
- property wl_eff¶
- property wl_range¶
- lightcurve_fitting.filters.extinction_law(freq, ebv, rv=3.1)[source]¶
A vectorized version of the Fitzpatrick (1999) extinction law from the
extinction
package- Parameters:
- freqarray-like
Frequencies in THz in the frame of the dust
- ebvarray-like
Selective extinction \(E(B-V)\) in magnitudes
- rvfloat, optional
Ratio of total to selective extinction \(R_V\). Default: 3.1.
- Returns:
- extinctionarray-like
Extinction factor \(10^{A/-2.5}\) at each input wavelength.
lightcurve_fitting.fitting module¶
- lightcurve_fitting.fitting.format_credible_interval(x, sigfigs=1, percentiles=(15.87, 50.0, 84.14), axis=0, varnames=None, units=None)[source]¶
Use LaTeX to format an equal-tailed credible interval with a given number of significant figures in the uncertainty
- Parameters:
- xarray-like
Data from which to calculate the credible interval
- sigfigsint, optional
Number of significant figures in the uncertainty. Default: 1
- percentilestuple, optional
Percentiles for the (lower, center, upper) of the credible interval. Default:
(15.87, 50., 84.14)
(median +/- 1σ)- axisint, optional
Axis of
x
along which to calculate the credible intervals. Default: 0- varnameslist, optional
Variable names to be equated with the credible intervals. Default: no variable names
- unitslist, optional
Units to be applied to the credible intervals. Default: no units
- Returns:
- paramtextsstr
The formatted credible intervals
- lightcurve_fitting.fitting.lightcurve_corner(lc, model, sampler_flatchain, model_kwargs=None, num_models_to_plot=100, lcaxis_posn=(0.7, 0.55, 0.2, 0.4), filter_spacing=1.0, tmin=None, tmax=None, t0_offset=None, save_plot_as='', ycol=None, textsize='medium', param_textsize='large', use_sigma=False, xscale='linear', filters_to_model=None, label_filters=True, lc_plot_kwargs=None, model_plot_kwargs=None)[source]¶
Plot the posterior distributions in a corner (pair) plot, with an inset showing the observed and model light curves.
- Parameters:
- lclightcurve_fitting.lightcurve.LC
Table of broadband photometry including columns “MJD”, “mag”, “dmag”, “filter”
- modellightcurve_fitting.models.Model
The model that was fit to the light curve.
- sampler_flatchainarray-like
2D array containing the aggregated MCMC chain histories
- model_kwargsdict, optional
DEPRECATED: Keyword arguments are now included in the model initialization
- num_models_to_plotint, optional
Number of model realizations to plot in the light curve inset. Default: 100
- lcaxis_posntuple, optional
Light curve inset position and size specification in figure units: (left, bottom, width, height)
- filter_spacingfloat, optional
Spacing between filters in the light curve inset, in units determined by the order of magnitude of the luminosities. Default: 1.
- tmin, tmaxfloat, optional
Starting and ending times for which to plot the models in the light curve inset. Default: determined by the time range of the observed light curve.
- t0_offsetfloat, optional
- Reference time for the explosion time in the corner plot and the horizontal axis of the light curve inset.
Default: the earliest explosion time in sampler_flatchain, rounded down.
- save_plot_asstr, optional
Filename to which to save the resulting plot
- ycolstr, optional
Quantity to plot on the light curve inset. Choices: “lum”, “flux”, or “absmag”. Default: model.output_quantity
- textsizestr, optional
Font size for the x- and y-axis labels, as well as the tick labels. Default: ‘medium’
- param_textsizestr, optional
Font size for the parameter text. Default: ‘large’
- use_sigmabool, optional
If True, treat the last parameter as an intrinsic scatter parameter that does not get passed to the model
- xscalestr, optional
Scale for the x-axis of the model plot. Choices: “linear” (default) or “log”.
- filters_to_modellist, set, optional
(Unique) list of filters for which to calculate the model light curves. Default: all filters in lc.
- label_filtersbool, optional
Print the name and offset for each filter to the right of the light curve axis. Default: True
- lc_plot_kwargsdict, optional
Dictionary of keyword arguments for the plotting of the observed light curve (passed to the lc.plot command)
- model_plot_kwargsdict, optional
Dictionary of keyword arguments for the plotting of the model light curves (passed to the plt.plot command)
- Returns:
- figmatplotlib.pyplot.Figure
Figure object containing the plot
- corner_axarray-like
Array of matplotlib.pyplot.Axes objects corresponding to the corner plot
- axmatplotlib.pyplot.Axes
Axes object for the light curve inset
- lightcurve_fitting.fitting.lightcurve_mcmc(lc, model, priors=None, p_min=None, p_max=None, p_lo=None, p_up=None, nwalkers=100, nsteps=1000, nsteps_burnin=1000, model_kwargs=None, show=False, save_plot_as='', save_sampler_as='', use_sigma=False, sigma_type='relative')[source]¶
Fit an analytical model to observed photometry using a Markov-chain Monte Carlo routine
- Parameters:
- lclightcurve_fitting.lightcurve.LC
Table of broadband photometry including columns “MJD”, “mag”, “dmag”, “filter”
- modellightcurve_fitting.models.Model
The model to fit to the light curve. Available models:
models.ShockCooling
,models.ShockCooling2
,models.ShockCooling3
,models.CompanionShocking
,models.CompanionShocking2
,models.CompanionShocking3
- priorslist, optional
Prior probability distributions for each model parameter. Available priors:
models.UniformPrior
(default),models.LogUniformPrior
,models.GaussianPrior
- p_min, p_maxlist, optional
DEPRECATED: Use priors instead
- p_lolist
Lower bounds on the starting guesses for each paramter
- p_uplist
Upper bounds on the starting guesses for each parameter
- nwalkersint, optional
Number of walkers (chains) for the MCMC routine. Default: 100
- nstepsint, optional
Number of steps (iterations) for the MCMC routine, excluding burn-in. Default: 1000
- nsteps_burninint, optional
Number of steps (iterations) for the MCMC routine during burn-in. Default: 1000
- model_kwargsdict, optional
DEPRECATED: Keyword arguments are now included in the model initialization
- showbool, optional
If True, plot and display the chain histories
- save_plot_asstr, optional
Save a plot of the chain histories to this filename
- save_sampler_asstr, optional
Save the aggregated chain histories to this filename
- use_sigmabool, optional
If True, treat the last parameter as an intrinsic scatter parameter that does not get passed to the model
- sigma_typestr, optional
If ‘relative’ (default), sigma will be in units of the individual photometric uncertainties. If ‘absolute’, sigma will be in units of the median photometric uncertainty.
- Returns:
- sampleremcee.EnsembleSampler
EnsembleSampler object containing the results of the fit
- lightcurve_fitting.fitting.lightcurve_model_plot(lc, model, sampler_flatchain, model_kwargs=None, num_models_to_plot=100, filter_spacing=1.0, tmin=None, tmax=None, ycol=None, textsize='medium', ax=None, mjd_offset=None, use_sigma=False, xscale='linear', filters_to_model=None, label_filters=True, lc_plot_kwargs=None, model_plot_kwargs=None)[source]¶
Plot the observed and model light curves.
- Parameters:
- lclightcurve_fitting.lightcurve.LC
Table of broadband photometry including columns “MJD”, “mag”, “dmag”, “filter”
- modellightcurve_fitting.models.Model
The model that was fit to the light curve.
- sampler_flatchainarray-like
2D array containing the aggregated MCMC chain histories
- model_kwargsdict, optional
DEPRECATED: Keyword arguments are now included in the model initialization.
- num_models_to_plotint, optional
Number of model realizations to plot in the light curve inset. Default: 100
- filter_spacingfloat, optional
Spacing between filters in the light curve inset, in units determined by the order of magnitude of the luminosities. Default: 1.
- tmin, tmaxfloat, optional
Starting and ending times for which to plot the models in the light curve inset. Default: determined by the time range of the observed light curve.
- ycolstr, optional
Quantity to plot on the light curve inset. Choices: “lum”, “flux”, or “absmag”. Default: model.output_quantity
- textsizestr, optional
Font size for the x- and y-axis labels, as well as the tick labels. Default: ‘medium’
- axmatplotlib.pyplot.Axes
Axis on which to plot the light curves
- mjd_offsetfloat, optional
Reference time on the horizontal axis of the light curve inset. Default: determined by the starting time of the model light curve.
- use_sigmabool, optional
If True, treat the last parameter as an intrinsic scatter parameter that does not get passed to the model
- xscalestr, optional
Scale for the x-axis. Choices: “linear” (default) or “log”.
- filters_to_modellist, set, optional
(Unique) list of filters for which to calculate the model light curves. Default: all filters in lc.
- label_filtersbool, optional
Print the name and offset for each filter to the right of the light curve axis. Default: True
- lc_plot_kwargsdict, optional
Dictionary of keyword arguments for the plotting of the observed light curve (passed to the lc.plot command)
- model_plot_kwargsdict, optional
Dictionary of keyword arguments for the plotting of the model light curves (passed to the plt.plot command)
lightcurve_fitting.lightcurve module¶
- class lightcurve_fitting.lightcurve.Arrow(hx, hy)[source]¶
Bases:
Path
An downward-pointing arrow-shaped
Path
, whose head has half-widthhx
and heighthy
(in units of length)
- class lightcurve_fitting.lightcurve.LC(*args, **kwargs)[source]¶
Bases:
Table
A broadband light curve, stored as an
astropy.table.Table
- Attributes:
- nondetSigmasfloat
Significance level implied by nondetections in the light curve. Default: 3σ
- groupbyset
Column names to group by when binning the light curve. Default:
{'filter', 'source'}
- markersdict
Mapping of some light curve property (default:
'source'
or'telescope'
) to marker shapes
- bin(delta=0.3, groupby=None)[source]¶
Bin the light curve by averaging points within
delta
days of each other- Parameters:
- deltafloat, optional
Bin size, in days. Default: 0.3 days
- groupbyset, optional
Column names to group by before binning. Default:
{'filter', 'source'}
- Returns:
- lclightcurve_fitting.lightcurve.LC
Binned light curve
- calcAbsMag(dm=None, extinction=None, hostext=None, ebv=None, rv=None, host_ebv=None, host_rv=None, redshift=None)[source]¶
Calculate the
'absmag'
column from the'mag'
column by correcting for distance and extinction- Parameters:
- dmfloat, optional
Distance modulus. Default: calculate from
redshift
.- extinctiondict, optional
Milky Way extinction coefficients \(A_λ\) for each filter. Default: calculate from
ebv
andrv
.- hostextdict, optional
Host galaxy extinction coefficients \(A_λ\) for each filter. Default: calculate from
host_ebv
andhost_rv
.- ebvfloat, optional
Milky Way selective extinction \(E(B-V)\), used if
extinction
is not given. Default: 0.- host_ebvfloat, optional
Host galaxy selective extinction \(E(B-V)\), used if
hostext
is not given. Default: 0.- rvfloat, optional
Ratio of total to selective Milky Way extinction \(R_V\), used with the
ebv
argument. Default: 3.1.- host_rvfloat, optional
Ratio of total to selective host-galaxy extinction \(R_V\), used with the
host_ebv
argument. Default: 3.1.- redshiftfloat, optional
Redshift of the host galaxy. Used to redshift the filters for host galaxy extinction. If no distance modulus is given, a redshift-dependent distance is calculated using the Planck18 cosmology. Default: 0.
- calcFlux(nondetSigmas=None, zp=None)[source]¶
Calculate the
'flux'
and'dflux'
columns from the'mag'
and'dmag'
columns- Parameters:
- nondetSigmasfloat, optional
Significance level implied by nondetections in the light curve. Default: 3σ
- zpfloat, array-like, optional
Array of zero points for each magnitude. Default: standard for each filter
- calcLum(nondetSigmas=None)[source]¶
Calculate the
'lum'
and'dlum'
columns from the'absmag'
and'dmag'
columns- Parameters:
- nondetSigmasfloat, optional
Significance level implied by nondetections in the light curve. Default: 3σ
- calcMag(nondetSigmas=None, zp=None)[source]¶
Calculate the
'mag'
and'dmag'
columns from the'flux'
and'dflux'
columns- Parameters:
- nondetSigmasfloat, optional
Significance level implied by nondetections in the light curve. Default: 3σ
- zpfloat, array-like, optional
Array of zero points for each magnitude. Default: standard for each filter
- calcPhase(rdsp=False, hours=False)[source]¶
Calculate the rest-frame
'phase'
column from'MJD'
,.meta['refmjd']
, and.meta['redshift']
- Parameters:
- rdspbool, optional
Define phase as rest-frame days since peak, rather than rest-frame days since explosion
- hoursbool, optional
Give the phase in rest-frame hours instead of rest-frame days
- findNondet(nondetSigmas=None)[source]¶
Add a boolean column
'nondet'
indicating flux measurements that are below the detection threshold- Parameters:
- nondetSigmasfloat, optional
Significance level implied by nondetections in the light curve. Default: 3σ
- findPeak(**criteria)[source]¶
Find the peak of the light curve and store it in
.meta['peakdate']
- Parameters:
- criteriadict, optional
Use only a subset of the light curve matching some criteria when calculating the peak date (stored in
.meta['peakcriteria']
)
- normalize_column_names()[source]¶
Rename any recognizable columns to their standard names for this package (see lightcurve.column_names).
- plot(xcol='phase', ycol='absmag', offset_factor=1.0, color='filter', marker=None, use_lines=False, normalize=False, fillmark=True, mjd_axis=True, appmag_axis=True, loc_mark=None, loc_filt=None, ncol_mark=1, lgd_filters=None, tight_layout=True, phase_hours=False, return_axes=False, frameon=True, **kwargs)[source]¶
Plot the light curve, with nondetections marked with a downward-pointing arrow
- Parameters:
- xcolstr, optional
Column to plot on the horizontal axis. Default:
'phase'
- ycolstr, optional
Column to plot on the vertical axis. Default:
'absmag'
- offset_factorfloat, optional
Increase or decrease the filter offsets by a constant factor. Default: 1.
- colorstr, optional
Column that controls the color of the lines and points. Default:
'filter'
- markerstr, optional
Column that controls the marker shape. Default:
'source'
or'telescope'
- use_linesbool, optional
Connect light curve points with lines. Default: False
- normalizebool, optional
Normalize all light curves to peak at 0. Default: False
- fillmarkbool, optional
Fill each marker with color. Default: True
- mjd_axisbool, optional
Plot MJD on the upper x-axis. Must have
.meta['redshift']
and.meta['refmjd']
. Default: True.- appmag_axisbool, optional
Plot extinction-corrected apparent magnitude on the right y-axis. Must have
.meta['dm']
. Default: True.- loc_mark, loc_filtstr, optional
Location for the marker and filter legends, respectively. Set to ‘none’ to omit them. Three new options are available: ‘above’, ‘above left’, and ‘above right’.
mjd_axis
and/orappmag_axis
must be used to add these legends. Otherwise runplt.legend()
after plotting for a single simple legend. Default: no legend.- ncol_markint, optional
Number of columns in the marker legend. Default: 1.
- lgd_filterslist, array-like, optional
Customize the arrangement of filters in the legend by providing a list of filters for each column.
None
can be used to leave a blank space in the column. Only filters given here will be used. The default arrangement shows all filters arranged by.system
(columns) and.offset
(rows).- tight_layoutbool, optional
Adjust the figure margins to look beautiful. Default: True.
- phase_hoursbool, optional
Plot the phase in units of rest-frame hours instead of rest-frame days. Default: False.
- return_axesbool, optional
Return the newly created axes if
mjd_axis=True
orappmag_axis=True
. Default: False.- kwargs
Keyword arguments matching column names in the light curve are used to specify a subset of points to plot. Additional keyword arguments passed to
matplotlib.pyplot.plot()
.
- Returns:
- topmatplotlib.pyplot.Axes, optional
The upper x-axis, if
mjd_axis=True
andreturn_axes=True
. Otherwise, None.- rightmatplotlib.pyplot.Axes, optional
The right y-axis, if
appmag_axis=True
andreturn_axes=True
. Otherwise, None.
- classmethod read(filepath, format='ascii', fill_values=None, **kwargs)[source]¶
Read and parse a data table and return as a Table.
This function provides the Table interface to the astropy unified I/O layer. This allows easily reading a file in many supported data formats using syntax such as:
>>> from astropy.table import Table >>> dat = Table.read('table.dat', format='ascii') >>> events = Table.read('events.fits', format='fits')
Get help on the available readers for
Table
using the``help()`` method:>>> Table.read.help() # Get help reading Table and list supported formats >>> Table.read.help('fits') # Get detailed help on Table FITS reader >>> Table.read.list_formats() # Print list of available formats
See also: https://docs.astropy.org/en/stable/io/unified.html
- Parameters:
- *argstuple, optional
Positional arguments passed through to data reader. If supplied the first argument is typically the input filename.
- formatstr
File format specifier.
- unitslist, dict, optional
List or dict of units to apply to columns
- descriptionslist, dict, optional
List or dict of descriptions to apply to columns
- **kwargsdict, optional
Keyword arguments passed through to data reader.
- Returns:
- out~astropy.table.Table
Table corresponding to file contents
Notes
The available built-in formats are:
Format
Read
Write
Auto-identify
Deprecated
ascii
Yes
Yes
No
ascii.aastex
Yes
Yes
No
ascii.basic
Yes
Yes
No
ascii.cds
Yes
No
No
ascii.commented_header
Yes
Yes
No
ascii.csv
Yes
Yes
Yes
ascii.daophot
Yes
No
No
ascii.ecsv
Yes
Yes
Yes
ascii.fast_basic
Yes
Yes
No
ascii.fast_commented_header
Yes
Yes
No
ascii.fast_csv
Yes
Yes
No
ascii.fast_no_header
Yes
Yes
No
ascii.fast_rdb
Yes
Yes
No
ascii.fast_tab
Yes
Yes
No
ascii.fixed_width
Yes
Yes
No
ascii.fixed_width_no_header
Yes
Yes
No
ascii.fixed_width_two_line
Yes
Yes
No
ascii.html
Yes
Yes
Yes
ascii.ipac
Yes
Yes
No
ascii.latex
Yes
Yes
Yes
ascii.mrt
Yes
Yes
No
ascii.no_header
Yes
Yes
No
ascii.qdp
Yes
Yes
Yes
ascii.rdb
Yes
Yes
Yes
ascii.rst
Yes
Yes
No
ascii.sextractor
Yes
No
No
ascii.tab
Yes
Yes
No
fits
Yes
Yes
Yes
hdf5
Yes
Yes
Yes
pandas.csv
Yes
Yes
No
pandas.fwf
Yes
No
No
pandas.html
Yes
Yes
No
pandas.json
Yes
Yes
No
parquet
Yes
Yes
Yes
votable
Yes
Yes
Yes
aastex
Yes
Yes
No
Yes
cds
Yes
No
No
Yes
csv
Yes
Yes
No
Yes
daophot
Yes
No
No
Yes
html
Yes
Yes
No
Yes
ipac
Yes
Yes
No
Yes
latex
Yes
Yes
No
Yes
mrt
Yes
Yes
No
Yes
rdb
Yes
Yes
No
Yes
Deprecated format names like
aastex
will be removed in a future version. Use the full name (e.g.ascii.aastex
) instead.
- where(**kwargs)[source]¶
Select the subset of a light curve matching some criteria, given as keyword arguments, e.g.,
colname=value
.- The keyword
colname
can be any of the following: a column in the table, in which case rows must match
value
in that columna column in the table +
_not
, in which case rows must not matchvalue
in that columna column in the table +
_min
, in which case rows must be >=value
in that columna column in the table +
_max
, in which case rows must be <=value
in that column
value
must match the data type of the columncolname
and can either be a single value or a list of values. Ifvalue
is a list, rows must match at least one of the values. Ifvalue
is a list andcolname
ends in_not
, rows must not match any of the values.- The keyword
- write(*args, **kwargs)[source]¶
Write this Table object out in the specified format.
This function provides the Table interface to the astropy unified I/O layer. This allows easily writing a file in many supported data formats using syntax such as:
>>> from astropy.table import Table >>> dat = Table([[1, 2], [3, 4]], names=('a', 'b')) >>> dat.write('table.dat', format='ascii')
Get help on the available writers for
Table
using the``help()`` method:>>> Table.write.help() # Get help writing Table and list supported formats >>> Table.write.help('fits') # Get detailed help on Table FITS writer >>> Table.write.list_formats() # Print list of available formats
The
serialize_method
argument is explained in the section on Table serialization methods.See also: https://docs.astropy.org/en/stable/io/unified.html
- Parameters:
- *argstuple, optional
Positional arguments passed through to data writer. If supplied the first argument is the output filename.
- formatstr
File format specifier.
- serialize_methodstr, dict, optional
Serialization method specifier for columns.
- **kwargsdict, optional
Keyword arguments passed through to data writer.
Notes
The available built-in formats are:
Format
Read
Write
Auto-identify
Deprecated
ascii
Yes
Yes
No
ascii.aastex
Yes
Yes
No
ascii.basic
Yes
Yes
No
ascii.commented_header
Yes
Yes
No
ascii.csv
Yes
Yes
Yes
ascii.ecsv
Yes
Yes
Yes
ascii.fast_basic
Yes
Yes
No
ascii.fast_commented_header
Yes
Yes
No
ascii.fast_csv
Yes
Yes
No
ascii.fast_no_header
Yes
Yes
No
ascii.fast_rdb
Yes
Yes
No
ascii.fast_tab
Yes
Yes
No
ascii.fixed_width
Yes
Yes
No
ascii.fixed_width_no_header
Yes
Yes
No
ascii.fixed_width_two_line
Yes
Yes
No
ascii.html
Yes
Yes
Yes
ascii.ipac
Yes
Yes
No
ascii.latex
Yes
Yes
Yes
ascii.mrt
Yes
Yes
No
ascii.no_header
Yes
Yes
No
ascii.qdp
Yes
Yes
Yes
ascii.rdb
Yes
Yes
Yes
ascii.rst
Yes
Yes
No
ascii.tab
Yes
Yes
No
fits
Yes
Yes
Yes
hdf5
Yes
Yes
Yes
jsviewer
No
Yes
No
pandas.csv
Yes
Yes
No
pandas.html
Yes
Yes
No
pandas.json
Yes
Yes
No
parquet
Yes
Yes
Yes
votable
Yes
Yes
Yes
votable.parquet
No
Yes
No
aastex
Yes
Yes
No
Yes
csv
Yes
Yes
No
Yes
html
Yes
Yes
No
Yes
ipac
Yes
Yes
No
Yes
latex
Yes
Yes
No
Yes
mrt
Yes
Yes
No
Yes
rdb
Yes
Yes
No
Yes
Deprecated format names like
aastex
will be removed in a future version. Use the full name (e.g.ascii.aastex
) instead.
- property zp¶
Returns an array of zero points for each filter in the
'filter'
column
- lightcurve_fitting.lightcurve.aux_axes(xfunc=None, yfunc=None, ax0=None, xfunc_args=None, yfunc_args=None)[source]¶
Add auxiliary axes to a plot that are linear transformations of the existing axes
- Parameters:
- xfuncfunction, optional
Function that transforms the lower x-axis to the upper x-axis. Default: do not add an upper x-axis.
- yfuncfunction, optional
Function that transforms the left y-axis to the right y-axis. Default: do not add a right y-axis.
- ax0matplotlib.pyplot.Axes, optional
Existing axes object. Default: use the current Axes.
- xfunc_args, yfunc_argsdict, optional
Keyword arguments for
xfunc
andyfunc
, respectively
- Returns:
- topmatplotlib.pyplot.Axes
The upper x-axis, if any. Otherwise, None.
- rightmatplotlib.pyplot.Axes
The right y-axis, if any. Otherwise, None.
- lightcurve_fitting.lightcurve.binflux(time, flux, dflux, delta=0.2, include_zero=True)[source]¶
Bin a light curve by averaging points within
delta
of each other in time- Parameters:
- time, flux, dfluxarray-like
Arrays of times, fluxes, and uncertainties comprising the observed light curve
- deltafloat, optional
Bin size, in the same units as
time
. Default: 0.2- include_zerobool, optional
Include data points with no error bar
- Returns:
- time, flux, dfluxarray-like
Binned arrays of times, fluxes, and uncertainties
- lightcurve_fitting.lightcurve.custom_legend(ax, handles, labels, top_axis=True, **kwargs)[source]¶
Add a legend to the axes with options for
loc='above'
,loc='above left'
, andloc='above right'
- Parameters:
- axmatplotlib.pyplot.Axes, matplotlib.pyplot.Figure
Axes or Figure object to which to add the legend
- handleslist of matplotlib.pyplot.Artist
A list of Artists (lines, patches) to be added to the legend
- labelslist of str
A list of labels to show next to the handles
- top_axisbool, optional
For legends above the top of the plot, add extra padding to the upper x-axis labels. Default: True.
- kwargs
Keyword arguments to be passed to
matplotlib.pyplot.legend()
- Returns:
- lgdmatplotlib.legend.Legend
The Legend object
- lightcurve_fitting.lightcurve.filter_legend(filts, offset_factor=1.0)[source]¶
Creates dummy artists and labels for the filter legend using the filter properties
- Parameters:
- filtsset, list, array-like
If a list or array of strings, the arrangement of filters in the legend, with columns in the first dimension and rows in the second. If a set of
Filter
objects, they will first be arranged usingfiltsetup()
.- offset_factorfloat, optional
Increase or decrease the filter offsets by a constant factor. Default: 1.
- Returns:
- lineslist of matplotlib.pyplot.Artist
A list of Artists (lines, patches) to be added to the legend
- labelslist of str
A list of labels to show next to the handles
- ncolint
Number of columns needed for the filter legend
- lightcurve_fitting.lightcurve.filtsetup(filts)[source]¶
Arrange filters in a grid according to their system (columns) and offset (rows)
- lightcurve_fitting.lightcurve.flux2mag(flux, dflux=array(nan), zp=0.0, nondet=None, nondetSigmas=3.0)[source]¶
Convert flux (and uncertainty) to magnitude (and uncertainty). Nondetections are converted to limiting magnitudes.
- Parameters:
- fluxfloat, array-like
Flux to be converted
- dfluxfloat, array-like, optional
Uncertainty on the flux to be converted. Default:
numpy.nan
- zpfloat, array-like, optional
Zero point to be applied to the magnitudes
- nondetarray-like
Boolean or index array indicating the nondetections among the fluxes. Default: no nondetections
- nondetSigmasfloat, optional
Significance level implied by nondetections in the light curve. Default: 3σ
- Returns:
- magfloat, array-like
Magnitude corresponding to the input flux
- dmagfloat, array-like
Uncertainty on the output magnitude
- lightcurve_fitting.lightcurve.mag2flux(mag, dmag=nan, zp=0.0, nondet=None, nondetSigmas=3.0)[source]¶
Convert magnitude (and uncertainty) to flux (and uncertainty). Nondetections are assumed to imply zero flux.
- Parameters:
- magfloat, array-like
Magnitude to be converted
- dmagfloat, array-like, optional
Uncertainty on the magnitude to be converted. Default:
numpy.nan
- zpfloat, array-like, optional
Zero point to be applied to the magnitudes
- nondetarray-like
Boolean or index array indicating the nondetections among the fluxes. Default: no nondetections
- nondetSigmasfloat, optional
Significance level implied by nondetections in the light curve. Default: 3σ
- Returns:
- fluxfloat, array-like
Flux corresponding to the input magnitude
- dfluxfloat, array-like
Uncertainty on the output flux
lightcurve_fitting.models module¶
- class lightcurve_fitting.models.BaseCompanionShocking(lc, redshift=0.0)[source]¶
Bases:
Model
The companion shocking model of Kasen (https://doi.org/10.1088/0004-637X/708/2/1025) plus the SiFTO SN Ia model.
As written by Hosseinzadeh et al. (https://doi.org/10.3847/2041-8213/aa8402), the shock component is defined by:
\(R_\mathrm{phot}(t) = (2700\,R_\odot) (M_c v_9^7)^{1/9} κ^{1/9} t^{7/9}\) (Eq. 1)
\(T_\mathrm{eff}(t) = (25\,\mathrm{kK}) a_{13}^{1/4} (M_c v_9^7)^{1/144} κ^{-35/144} t^{37/72}\) (Eq. 2)
The SiFTO model (https://doi.org/10.1086/588518) is currently only available in the UBVgri filters.
- Parameters:
- lclightcurve_fitting.lightcurve.LC, optional
The light curve to which the model will be fit. Only used for lc.meta[‘redshift’] if z is not given.
- redshiftfloat, optional
The redshift between blackbody source and the observed filters. Default: 0.
- Attributes:
- zfloat
The redshift between blackbody source and the observed filters
- siftoastropy.table.Table
A copy of the SiFTO model scaled to match the observed peak luminosity in each filter
- companion_shocking(t_in, f, t_exp, a13, Mc_v9_7, kappa=1.0)[source]¶
Evaluate the shock component only at a range of times and filters
- Parameters:
- t_infloat, array-like
Time in days
- flightcurve_fitting.filter.Filter, array-like
Filters for which to calculate the model
- t_expfloat, array-like
The explosion epoch
- a13float, array-like
The binary separation in \(10^{13}\) cm
- Mc_v9_7float, array-like
The product \(M_c v_9^7\), where \(M_c\) is the ejecta mass in Chandrasekhar masses and \(v_9\) is the ejecta velocity in units of \(10^9\) cm/s
- kappafloat, array-like
The ejecta opacity in units of the electron scattering opacity (0.34 cm^2/g)
- Returns:
- y_fitarray-like
The filtered model light curves
- stretched_sifto(t_in, f, t_peak, stretch, dtU=None, dti=None)[source]¶
The SiFTO SN Ia model (https://doi.org/10.1086/588518), offset and stretched by the input parameters
The SiFTO model is currently only available in the UBVgri filters. We assume unfiltered photometry can be modeled as r.
- Parameters:
- t_infloat, array-like
Time in days
- flightcurve_fitting.filter.Filter, array-like
Filters for which to calculate the model
- t_peakfloat, array-like
The epoch of maximum light for the SiFTO model
- stretchfloat, array-like
The stretch for the SiFTO model
- dtU, dtifloat, array-like
Time offsets for the U- and i-band models relative to the other bands
- Returns:
- y_fitarray-like
The filtered model light curves
- static t_max(p)[source]¶
The maximum time at which the model is valid
This is the last epoch at which the stretched SiFTO model is computed
- static t_min(p)[source]¶
The minimum time at which the model is valid
This is the first epoch at which the stretched SiFTO model is computed
- static temperature_radius(t_in, t_exp, a13, Mc_v9_7, kappa=1.0)[source]¶
Evaluate the color temperature and photospheric radius of the shock component as a function of time
- Parameters:
- t_infloat, array-like
Time in days
- t_expfloat, array-like
The explosion epoch
- a13float, array-like
The binary separation in \(10^{13}\) cm
- Mc_v9_7float, array-like
The product \(M_c v_9^7\), where \(M_c\) is the ejecta mass in Chandrasekhar masses and \(v_9\) is the ejecta velocity in units of \(10^9\) cm/s
- kappafloat, array-like
The ejecta opacity in units of the electron scattering opacity (0.34 cm^2/g)
- Returns:
- T_kasenarray-like
The model blackbody temperatures in units of kilokelvins
- R_kasenarray-like
The model blackbody radii in units of 1000 solar radii
- class lightcurve_fitting.models.BaseShockCooling(lc=None, redshift=0.0, n=1.5, RW=False)[source]¶
Bases:
Model
The shock cooling model of Sapir & Waxman (https://doi.org/10.3847/1538-4357/aa64df).
\(T(t) = \frac{T_\mathrm{col}}{T_\mathrm{ph}} T_0 \left(\frac{v_s^2 t^2}{f_ρ M κ}\right)^{ε_1} \frac{R^{1/4}}{κ^{1/4}} t^{-1/2}\) (Eq. 23)
\(L(t) = A \exp\left[-\left(\frac{a t}{t_\mathrm{tr}}\right)^α\right] L_0 \left(\frac{v_s t^2}{f_ρ M κ}\right)^{-ε_2} \frac{v_s^2 R}{κ}\) (Eq. 18-19)
\(t_\mathrm{tr} = (19.5\,\mathrm{d}) \left(\frac{κ * M_\mathrm{env}}{v_s} \right)^{1/2}\) (Eq. 20)
- Parameters:
- lclightcurve_fitting.lightcurve.LC, optional
The light curve to which the model will be fit. Only used to get the redshift if redshift is not given.
- redshiftfloat, optional
The redshift between blackbody source and the observed filters. Default: 0.
- nfloat, optional
The polytropic index of the progenitor. Must be either 1.5 (default) or 3.
- RWbool, optional
Reduce the model to the simpler form of Rabinak & Waxman (https://doi.org/10.1088/0004-637X/728/1/63). Default: False.
- Attributes:
- zfloat
The redshift between blackbody source and the observed filters
- nfloat
The polytropic index of the progenitor
- Afloat
Coefficient on the luminosity suppression factor (Eq. 19)
- afloat
Coefficient on the transparency timescale (Eq. 19)
- alphafloat
Exponent on the transparency timescale (Eq. 19)
- epsilon_1float
Exponent in the temperature expression (Eq. 23)
- epsilon_2float
Exponent in the luminosity expression (Eq. 18)
- L_0float
Coefficient on the luminosity expression in erg/s (Eq. 18)
- T_0float
Coefficient on the temperature expression in eV (Eq. 23)
- Tph_to_Tcolfloat
Ratio of the color temperature to the photospheric temperature
- epsilon_Tfloat
Exponent on time in the temperature expression
- epsilon_Lfloat
Exponent on time in the luminosity expression
- RWbool
Reduce the model to the simpler form of Rabinak & Waxman (https://doi.org/10.1088/0004-637X/728/1/63)
- static t_max(p, kappa=1.0)[source]¶
The maximum time at which the model is valid
\(t_\mathrm{max} = (7.4\,\mathrm{d}) \left(\frac{R}{κ}\right)^{0.55} + t_\mathrm{exp}\) (Eq. 24)
- static t_min(p, kappa=1.0)[source]¶
The minimum time at which the model is valid
\(t_\mathrm{min} = (0.2\,\mathrm{d}) \frac{R}{v_s} \max\left[0.5, \frac{R^{0.4}}{(f_ρ M κ)^{0.2} v_s^{-0.7}}\right] + t_\mathrm{exp}\) (Eq. 17)
- temperature_radius(t_in, v_s, M_env, f_rho_M, R, t_exp=0.0, kappa=1.0)[source]¶
Evaluate the color temperature and photospheric radius as a function of time for a set of parameters
- Parameters:
- t_infloat, array-like
Time in days
- v_sfloat, array-like
The shock speed in \(10^{8.5}\) cm/s
- M_envfloat, array-like
The envelope mass in solar masses
- f_rho_Mfloat, array-like
The product \(f_ρ M\), where \(f_ρ\) is a numerical factor of order unity that depends on the inner envelope structure and \(M\) is the ejecta mass in solar masses
- Rfloat, array-like
The progenitor radius in \(10^{13}\) cm
- t_expfloat, array-like
The explosion epoch
- kappafloat, array-like
The ejecta opacity in units of the electron scattering opacity (0.34 cm^2/g)
- Returns:
- T_Karray-like
The model blackbody temperatures in units of kilokelvins
- R_bbarray-like
The model blackbody radii in units of 1000 solar radii
- class lightcurve_fitting.models.CompanionShocking(lc, redshift=0.0)[source]¶
Bases:
BaseCompanionShocking
The companion shocking model of Kasen (https://doi.org/10.1088/0004-637X/708/2/1025) plus the SiFTO SN Ia model.
This version of the model includes factors on the r and i SiFTO models, and a factor on the U shock component.
- evaluate(t_in, f, t_exp, a13, Mc_v9_7, t_peak, stretch, rr=1.0, ri=1.0, rU=1.0, kappa=1.0)[source]¶
Evaluate this model at a range of times and filters
- Parameters:
- t_infloat, array-like
Time in days
- flightcurve_fitting.filter.Filter, array-like
Filters for which to calculate the model
- t_expfloat, array-like
The explosion epoch
- a13float, array-like
The binary separation in \(10^{13}\) cm
- Mc_v9_7float, array-like
The product \(M_c v_9^7\), where \(M_c\) is the ejecta mass in Chandrasekhar masses and \(v_9\) is the ejecta velocity in units of \(10^9\) cm/s
- t_peakfloat, array-like
The epoch of maximum light for the SiFTO model
- stretchfloat, array-like
The stretch for the SiFTO model
- rrfloat, array-like
A scale factor for the r-band SiFTO model
- rifloat, array-like
A scale factor for the i-band SiFTO model
- rUfloat, aray-like
A scale factor for the U-band of the shock component
- kappafloat, array-like
The ejecta opacity in units of the electron scattering opacity (0.34 cm^2/g)
- Returns:
- y_fitarray-like
The filtered model light curves
- input_names = ['t_0', 'a', 'M v^7', 't_\\mathrm{max}', 's', 'r_r', 'r_i', 'r_U']¶
A list of the parameter names
- units = [Unit("d"), <Quantity 1.e+13 cm>, <Quantity 1.e+63 cm7 M_chandra / s7>, Unit("d"), Unit(dimensionless), Unit(dimensionless), Unit(dimensionless), Unit(dimensionless)]¶
A list of the units for each parameter
- class lightcurve_fitting.models.CompanionShocking2(lc, redshift=0.0)[source]¶
Bases:
BaseCompanionShocking
The companion shocking model of Kasen (https://doi.org/10.1088/0004-637X/708/2/1025) plus the SiFTO SN Ia model.
This version of the model includes time offsets for the U and i SiFTO models.
- evaluate(t_in, f, t_exp, a13, Mc_v9_7, t_peak, stretch, dtU=0.0, dti=0.0, kappa=1.0)[source]¶
Evaluate this model at a range of times and filters
- Parameters:
- t_infloat, array-like
Time in days
- flightcurve_fitting.filter.Filter, array-like
Filters for which to calculate the model
- t_expfloat, array-like
The explosion epoch
- a13float, array-like
The binary separation in \(10^{13}\) cm
- Mc_v9_7float, array-like
The product \(M_c v_9^7\), where \(M_c\) is the ejecta mass in Chandrasekhar masses and \(v_9\) is the ejecta velocity in units of \(10^9\) cm/s
- t_peakfloat, array-like
The epoch of maximum light for the SiFTO model
- stretchfloat, array-like
The stretch for the SiFTO model
- dtU, dtifloat, array-like
Time offsets for the U- and i-band SiFTO models relative to the other bands
- kappafloat, array-like
The ejecta opacity in units of the electron scattering opacity (0.34 cm^2/g)
- Returns:
- y_fitarray-like
The filtered model light curves
- input_names = ['t_0', 'a', 'M v^7', 't_\\mathrm{max}', 's', '\\Delta t_U', '\\Delta t_i']¶
A list of the parameter names
- units = [Unit("d"), <Quantity 1.e+13 cm>, <Quantity 1.e+63 cm7 M_chandra / s7>, Unit("d"), Unit(dimensionless), Unit("d"), Unit("d")]¶
A list of the units for each parameter
- class lightcurve_fitting.models.CompanionShocking3(lc, redshift=0.0)[source]¶
Bases:
BaseCompanionShocking
The companion shocking model of Kasen (https://doi.org/10.1088/0004-637X/708/2/1025) plus the SiFTO SN Ia model.
This version of the model includes time offsets for the U and i SiFTO models, as well as the viewing angle dependence parametrized by Brown et al. (https://doi.org/10.1088/0004-637X/749/1/18).
- evaluate(t_in, f, t_exp, a13, theta, t_peak, stretch, dtU, dti, kappa=1.0)[source]¶
Evaluate this model at a range of times and filters
- Parameters:
- t_infloat, array-like
Time in days
- flightcurve_fitting.filter.Filter, array-like
Filters for which to calculate the model
- t_expfloat, array-like
The explosion epoch
- a13float, array-like
The binary separation in \(10^{13}\) cm
- thetafloat, array-like
The angle between the binary axis and the line to the observer in degrees. 0° means the binary companion is along the line of sight.
- t_peakfloat, array-like
The epoch of maximum light for the SiFTO model
- stretchfloat, array-like
The stretch for the SiFTO model
- dtU, dtifloat, array-like
Time offsets for the U- and i-band SiFTO models relative to the other bands
- kappafloat, array-like
The ejecta opacity in units of the electron scattering opacity (0.34 cm^2/g)
- Returns:
- y_fitarray-like
The filtered model light curves
- input_names = ['t_0', 'a', '\\theta', 't_\\mathrm{max}', 's', '\\Delta t_U', '\\Delta t_i']¶
A list of the parameter names
- units = [Unit("d"), <Quantity 1.e+13 cm>, Unit("deg"), Unit("d"), Unit(dimensionless), Unit("d"), Unit("d")]¶
A list of the units for each parameter
- class lightcurve_fitting.models.GaussianPrior(p_min=-inf, p_max=inf, mean=0.0, stddev=1.0)[source]¶
Bases:
Prior
A Gaussian prior centered at mean with standard deviation stddev, i.e., \(\frac{dP}{dp} \propto \exp \left( \frac{(p - \mu)^2}{2 \sigma^2} \right)\)
- class lightcurve_fitting.models.LogUniformPrior(p_min=0.0, p_max=inf)[source]¶
Bases:
Prior
A uniform prior in the logarithm of the parameter, i.e., \(\frac{dP}{dp} \propto \frac{1}{p}\)
- class lightcurve_fitting.models.Model(lc=None, redshift=0.0)[source]¶
Bases:
object
An analytical model, defined by a function and its parameters
- property axis_labels¶
Axis labels for each paramter (including name and unit)
- input_names = []¶
A list of the parameter names
- log_likelihood(lc, p, use_sigma=False, sigma_type='relative')[source]¶
The log of the likelihood function of the model given data contained in lc and parameters contained in p
- Parameters:
- lclightcurve_fitting.lightcurve.LC
Table of broadband photometry including columns “MJD”, “mag”, “dmag”, “filter”
- parray-like
Model parameters for which to calculate the likelihood. The first dimension should have the same length as the number of parameters, or one more if sigma is used.
- use_sigmabool, optional
If True, treat the last parameter as an intrinsic scatter parameter that does not get passed to the model
- sigma_typestr, optional
If ‘relative’ (default), sigma will be in units of the individual photometric uncertainties. If ‘absolute’, sigma will be in units of the median photometric uncertainty.
- Returns:
- log_likelihoodfloat, array-like
The natural log of the likelihood function. If p is 1D, log_likelihood is a float. Otherwise, its shape is determined by the remaining dimensions of p.
- property nparams¶
The number of parameters in the model
- output_quantity = 'lum'¶
Quantity output by the model: ‘lum’ or ‘flux’
- units = []¶
A list of the units for each parameter
- class lightcurve_fitting.models.ShockCooling(lc=None, redshift=0.0, n=1.5, RW=False)[source]¶
Bases:
BaseShockCooling
The shock cooling model of Sapir & Waxman (https://doi.org/10.3847/1538-4357/aa64df).
This version of the model is written in terms of physical parameters \(v_s, M_\mathrm{env}, f_ρ M, R\).
- evaluate(t_in, f, v_s, M_env, f_rho_M, R, t_exp=0.0, kappa=1.0)[source]¶
Evaluate this model at a range of times and filters
- Parameters:
- t_infloat, array-like
Time in days
- flightcurve_fitting.filter.Filter, array-like
Filters for which to calculate the model
- v_sfloat, array-like
The shock speed in \(10^{8.5}\) cm/s
- M_envfloat, array-like
The envelope mass in solar masses
- f_rho_Mfloat, array-like
The product \(f_ρ M\), where “\(f_ρ\) is a numerical factor of order unity that depends on the inner envelope structure” and \(M\) is the ejecta mass in solar masses
- Rfloat, array-like
The progenitor radius in \(10^{13}\) cm
- t_expfloat, array-like, optional
The explosion epoch. Default: 0.
- kappafloat, array-like, optional
The ejecta opacity in units of the electron scattering opacity (0.34 cm^2/g). Default: 1.
- Returns:
- y_fitarray-like
The filtered model light curves
- input_names = ['v_\\mathrm{s*}', 'M_\\mathrm{env}', 'f_\\rho M', 'R', 't_0']¶
A list of the parameter names
- units = [<Quantity 3.16227766e+08 cm / s>, Unit("solMass"), Unit("solMass"), <Quantity 1.e+13 cm>, Unit("d")]¶
A list of the units for each parameter
- class lightcurve_fitting.models.ShockCooling2(lc=None, redshift=0.0, n=1.5, RW=False)[source]¶
Bases:
BaseShockCooling
The shock cooling model of Sapir & Waxman (https://doi.org/10.3847/1538-4357/aa64df).
This version of the model is written in terms of scaling parameters \(T_1, L_1, t_\mathrm{tr}\):
\(T(t) = T_1 t^{ε_T}\) and \(L(t) = L_1 t^{ε_L} \exp\left[-\left(\frac{a t}{t_\mathrm{tr}}\right)^α\right]\)
- evaluate(t_in, f, T_1, L_1, t_tr, t_exp=0.0)[source]¶
Evaluate this model at a range of times and filters
- Parameters:
- t_infloat, array-like
Time in days
- flightcurve_fitting.filter.Filter, array-like
Filters for which to calculate the model
- T_1float, array-like
The blackbody temperature 1 day after explosion in kilokelvins
- L_1float, array-like
The approximate blackbody luminosity 1 day after explosion in \(10^{42}\) erg/s
- t_trfloat, array-like
The time at which the envelope becomes transparent in rest-frame days
- t_expfloat, array-like
The explosion epoch
- Returns:
- y_fitarray-like
The filtered model light curves
- input_names = ['T_1', 'L_1', 't_\\mathrm{tr}', 't_0']¶
A list of the parameter names
- t_max(p, kappa=1.0)[source]¶
The maximum time at which the model is valid
\(t_\mathrm{max} = \left(\frac{8.12\,\mathrm{kK}}{T_1}\right)^{1/(2 ε_1 - 0.5)} + t_\mathrm{exp}\)
- static t_min(p, kappa=1.0)[source]¶
The minimum time at which the model is valid
This expression cannot be translated to the parameters of
ShockCooling2
- units = [Unit("kK"), <Quantity 1.e+42 erg / s>, Unit("d"), Unit("d")]¶
A list of the units for each parameter
- class lightcurve_fitting.models.ShockCooling3(lc=None, redshift=0.0, n=1.5, RW=False)[source]¶
Bases:
BaseShockCooling
The shock cooling model of Sapir & Waxman (https://doi.org/10.3847/1538-4357/aa64df).
This version of the model is written in terms of physical parameters \(v_s, M_\mathrm{env}, f_ρ M, R\) and includes distance \(d_L\) and reddening \(E(B-V)\) as free parameters.
- evaluate(t_in, f, v_s, M_env, f_rho_M, R, dist, ebv=0.0, t_exp=0.0, kappa=1.0)[source]¶
Evaluate this model at a range of times and filters
- Parameters:
- t_infloat, array-like
Time in days
- flightcurve_fitting.filter.Filter, array-like
Filters for which to calculate the model
- v_sfloat, array-like
The shock speed in \(10^{8.5}\) cm/s
- M_envfloat, array-like
The envelope mass in solar masses
- f_rho_Mfloat, array-like
The product \(f_ρ M\), where “\(f_ρ\) is a numerical factor of order unity that depends on the inner envelope structure” and \(M\) is the ejecta mass in solar masses
- Rfloat, array-like
The progenitor radius in \(10^{13}\) cm
- distfloat, array-like
The luminosity distance in Mpc
- ebvfloat, array-like, optional
The reddening \(E(B-V)\) to apply to the blackbody spectrum before integration. Default: 0.
- t_expfloat, array-like, optional
The explosion epoch. Default: 0.
- kappafloat, array-like, optional
The ejecta opacity in units of the electron scattering opacity (0.34 cm^2/g). Default: 1.
- Returns:
- y_fitarray-like
The filtered model light curves
- input_names = ['v_\\mathrm{s*}', 'M_\\mathrm{env}', 'f_\\rho M', 'R', 'd_L', 'E(B-V)', 't_0']¶
A list of the parameter names
- output_quantity = 'flux'¶
Quantity output by the model: ‘lum’ or ‘flux’
- static t_max(p, kappa=1.0)[source]¶
The maximum time at which the model is valid
\(t_\mathrm{max} = (7.4\,\mathrm{d}) \left(\frac{R}{κ}\right)^{0.55} + t_\mathrm{exp}\) (Eq. 24)
- static t_min(p, kappa=1.0)[source]¶
The minimum time at which the model is valid
\(t_\mathrm{min} = (0.2\,\mathrm{d}) \frac{R}{v_s} \max\left[0.5, \frac{R^{0.4}}{(f_ρ M κ)^{0.2} v_s^{-0.7}}\right] + t_\mathrm{exp}\) (Eq. 17)
- units = [<Quantity 3.16227766e+08 cm / s>, Unit("solMass"), Unit("solMass"), <Quantity 1.e+13 cm>, Unit("Mpc"), Unit("mag"), Unit("d")]¶
A list of the units for each parameter
- class lightcurve_fitting.models.ShockCooling4(lc=None, redshift=0.0)[source]¶
Bases:
Model
The shock cooling model of Morag, Sapir, & Waxman (https://doi.org/10.1093/mnras/stad899).
\(L(\tilde{t}) = L_\mathrm{br}\left\{\tilde{t}^{-4/3} + 0.9\exp\left[-\left(\frac{2.0t}{t_\mathrm{tr}}\right)^{0.5}\right] \tilde{t}^{-0.17}\right\}\) (Eq. A1)
\(T_\mathrm{col}(\tilde{t}) = T_\mathrm{col,br} \min(0.97\tilde{t}^{-1/3}, \tilde{t}^{-0.45})\) (Eq. A2)
\(\tilde{t} \equiv \frac{t}{t_\mathrm{br}}\), where \(t_\mathrm{br} = (0.86\,\mathrm{h}) R^{1.26} v_\mathrm{s*}^{-1.13} (f_\rho M \kappa)^{-0.13}\) (Eq. A5)
\(L_\mathrm{br} = (3.69 \times 10^{42}\,\mathrm{erg}\,\mathrm{s}^{-1}) R^{0.78} v_\mathrm{s*}^{2.11} (f_\rho M)^{0.11} \kappa^{-0.89}\) (Eq. A6)
\(T_\mathrm{col,br} = (8.19\,\mathrm{eV}) R^{-0.32} v_\mathrm{s*}^{0.58} (f_\rho M)^{0.03} \kappa^{-0.22}\) (Eq. A7)
\(t_\mathrm{tr} = (19.5\,\mathrm{d}) \sqrt{\frac{\kappa M}{v_\mathrm{s*}}}\) (Eq. A9)
- Parameters:
- lclightcurve_fitting.lightcurve.LC, optional
The light curve to which the model will be fit. Only used to get the redshift if redshift is not given.
- redshiftfloat, optional
The redshift between blackbody source and the observed filters. Default: 0.
- Attributes:
- zfloat
The redshift between blackbody source and the observed filters
- nfloat
The polytropic index of the progenitor
- Afloat
Coefficient on the luminosity suppression factor (Eq. A1)
- afloat
Coefficient on the transparency timescale (Eq. A1)
- alphafloat
Exponent on the transparency timescale (Eq. A1)
- L_br_0float
Coefficient on the luminosity expression in erg/s (Eq. A6)
- T_col_br_0float
Coefficient on the temperature expression in eV (Eq. A7)
- t_min_0float
Coefficient on the minimum validity time in days (Eq. A3)
- t_br_0float
Coefficient on the \(\tilde{t}\) timescale in days (Eq. A5)
- t_07eV_0float
Coefficient on the time at which the ejecta reach 0.7 eV in days (Eq. A8)
- t_tr_0float
Coefficient on the transparency timescale in days (Eq. A9)
- evaluate(t_in, f, v_s, M_env, f_rho_M, R, t_exp=0.0, kappa=1.0)[source]¶
Evaluate this model at a range of times and filters
- Parameters:
- t_infloat, array-like
Time in days
- flightcurve_fitting.filter.Filter, array-like
Filters for which to calculate the model
- v_sfloat, array-like
The shock speed in \(10^{8.5}\) cm/s
- M_envfloat, array-like
The envelope mass in solar masses
- f_rho_Mfloat, array-like
The product \(f_ρ M\), where “\(f_ρ\) is a numerical factor of order unity that depends on the inner envelope structure” and \(M\) is the ejecta mass in solar masses
- Rfloat, array-like
The progenitor radius in \(10^{13}\) cm
- t_expfloat, array-like, optional
The explosion epoch. Default: 0.
- kappafloat, array-like, optional
The ejecta opacity in units of the electron scattering opacity (0.34 cm^2/g). Default: 1.
- Returns:
- y_fitarray-like
The filtered model light curves
- input_names = ['v_\\mathrm{s*}', 'M_\\mathrm{env}', 'f_\\rho M', 'R', 't_0']¶
A list of the parameter names
- t_max(p, kappa=1.0)[source]¶
The maximum time at which the model is valid
\(t_\mathrm{max} = \min(t_\mathrm{0.7\,eV}, 0.5 t_\mathrm{tr})\) (Eq. A3)
\(t_\mathrm{0.7\,eV} = (6.86\,\mathrm{d}) R^{0.56} v_\mathrm{s*}^{0.16} \kappa^{-0.61} (f_\rho M)^{-0.06}\) (Eq. A8)
\(t_\mathrm{tr} = (19.5\,\mathrm{d}) \sqrt{\frac{\kappa M}{v_\mathrm{s*}}}\) (Eq. A9)
- t_min(p, kappa=1.0)[source]¶
The minimum time at which the model is valid
\(t_\mathrm{min} = (17\,\mathrm{min}) R + t_\mathrm{exp}\) (Eq. A3)
- units = [<Quantity 3.16227766e+08 cm / s>, Unit("solMass"), Unit("solMass"), <Quantity 1.e+13 cm>, Unit("d")]¶
A list of the units for each parameter
- class lightcurve_fitting.models.UniformPrior(p_min=-inf, p_max=inf)[source]¶
Bases:
Prior
A uniform prior in the parameter, i.e., \(\frac{dP}{dp} \propto 1\)
- lightcurve_fitting.models.blackbody_to_filters(filters, T, R, z=0.0, cutoff_freq=inf, ebv=0.0)[source]¶
The average spectral luminosity density (\(L_ν\)) of a blackbody as observed through one or more filters
- Parameters:
- filterslightcurve_fitting.filters.Filter, array-like
One or more broadband filters
- Tfloat, array-like
Blackbody temperatures in kilokelvins
- Rfloat, array-like
Blackbody radii in units of 1000 solar radii
- zfloat
Redshift between the blackbody source and the observed filters
- cutoff_freqfloat, optional
Cutoff frequency (in terahertz) for a modified blackbody spectrum (see https://doi.org/10.3847/1538-4357/aa9334)
- ebvfloat, array-like, optional
Selective extinction E(B-V) in magnitudes, evaluated using a Fitzpatrick (1999) extinction law with R_V=3.1. Its shape must be broadcastable to T and R. Default: 0.
- Returns:
- y_fitfloat, array-like
The average spectral luminosity density in each filter in watts per hertz
- lightcurve_fitting.models.format_unit(unit)[source]¶
Use LaTeX to format a physical unit, which may consist of an order of magnitude times a base unit
- Parameters:
- unitastropy.units.Unit, astropy.units.Quantity
The unit to format
- Returns:
- unit_strstr
Formatted unit
- lightcurve_fitting.models.planck(nu, T, R, dT=0.0, dR=0.0, cov=0.0)[source]¶
The Planck spectrum for a blackbody source, including propagation of uncertainties
- Parameters:
- nufloat, array-like
Frequency in terahertz
- Tfloat, array-like
Temperature in kilokelvins
- Rfloat, array-like
Radius in units of 1000 solar radii
- dTfloat, array-like, optional
Uncertainty in the temperature in kilokelvins
- dRfloat, array-like, optional
Uncertainty in the radius in units of 1000 solar radii
- covfloat, array-like, optional
The covariance between the temperature and radius
- Returns:
- Lnufloat, array-like
The spectral luminosity density (\(L_ν\)) of the source in watts per hertz
- dLnufloat, array-like
The uncertainty in the spectral luminosity density in watts per hertz
- lightcurve_fitting.models.planck_fast(nu, T, R, cutoff_freq=inf)[source]¶
The Planck spectrum for a blackbody source
\(L_ν = 4 π R^2 \frac{2 π h ν^3 c^{-2}}{\exp(h ν / k_B T) - 1}\)
- Parameters:
- nufloat, array-like
Frequency in terahertz
- Tfloat, array-like
Temperature in kilokelvins
- Rfloat, array-like
Radius in units of 1000 solar radii
- cutoff_freqfloat, optional
Cutoff frequency (in terahertz) for a modified blackbody spectrum (see https://doi.org/10.3847/1538-4357/aa9334)
- Returns:
- float, array-like
The spectral luminosity density (\(L_ν\)) of the source in watts per hertz
lightcurve_fitting.speccal module¶
- lightcurve_fitting.speccal.calibrate_spectra(spectra, lc, filters=None, order=0, subtract_percentile=None, max_extrapolate=1.0, show=False)[source]¶
Calibrate a set of spectra to an observed broadband light curve.
Each calibrated spectrum is saved to a text file that corresponds to the original filename prefixed by
photcal_
.- Parameters:
- spectralist
List of filenames containing the spectra
- lclightcurve_fitting.lightcurve.LC
Photometry table containing the observed light curve
- filterslist, optional
Only use this subset of filters for calibration
- orderint, optional
Polynomial order for the calibration function. Default: 0 (constant factor)
- subtract_percentilefloat, optional
Subtract flux corresponding to this percentile of the spectrum before calibration. Default: no subtraction
- max_extrapolatefloat, optional
Assume constant flux in a filter for this many days after the last observed point. Default: 1 day.
- showbool, optional
Plot the observed light curve and the uncalibrated and calibrated spectra, and ask whether to save the results
- Returns:
- figmatplotlib.pyplot.Figure
If show=True, the Figure object is returned
- lightcurve_fitting.speccal.convert_spectrum_units(wl, flux, hdr, default_bunit='erg / (Angstrom cm2 s)', default_cunit='Angstrom')[source]¶
Convert a spectrum to standard units, if information is available in the header to establish units (BUNIT & CUNIT1)
- Parameters:
- wl: array-like
Input wavelengths
- flux: array-like
Input fluxes
- hdr: dict-like
Metadata from which to determine the units of wl and flux
- default_bunit: str, optional
Units of the output flux. Default: ‘erg / (Angstrom cm2 s)’
- default_cunit: str, optional
Units of the output wavelengths. Default: ‘Angstrom’
- Returns:
- wl: numpy.ndarray
Input wavelengths converted to default_cunit
- flux: numpy.ndarray
Input fluxes converted to default_bunit
- lightcurve_fitting.speccal.create_wiserep_tsv(specpaths, wiserep_dir, verbose=False, instruments=None, date_fmt='iso')[source]¶
Prepares a TSV file for uploading spectra to WISeREP (see https://www.wiserep.org/content/wiserep-getting-started).
Also copies all the spectrum files to a single directory, and converts any FITS files to ASCII.
- Parameters:
- specpaths: list
List of paths to the spectrum files. Optionally, data quality can be specified in a tuple, e.g., [(filename1.fits, 3), (filename2.txt, 2), …], where the second element is 1 (low), 2 (medium), or 3 (high).
- wiserep_dir: str
Directory in which to collect the spectra. If it already exists, it will be deleted and remade. The TSV file will be saved alongside this directory using the same name but with a ‘.tsv’ extension.
- verbose: bool, optional
If True, print a message when any file is copied or created.
- instruments: dict, optional
A dictionary of known instruments. The keys are the instrument from the header and the values are the wiserep instrument ids from https://www.wiserep.org/aux. By default the script builds this dictionary on the fly.
- date_fmt: str, optional
Format for the Obs-date column in the output. WISeREP accepts ‘iso’ (default) or ‘jd’.
- lightcurve_fitting.speccal.readOSCspec(filepath)[source]¶
Read spectra from a JSON file in the Open Astronomy Catalog schema (https://github.com/astrocatalogs/schema)
- Parameters:
- filepathstr
Path to the JSON file containing the spectra
- Returns:
- filenameslist
List of filenames corresponding to each spectrum
- timeslist
List of
astropy.time.Time
when each spectrum was observed- tellist
List of telescope names (if given) where each spectrum was observed
- instlist
List of instrument names (if given) with which each spectrum was observed
- wllist
List of wavelength arrays for each spectrum
- fxlist
List of flux arrays for each spectrum
- scaleslist
List of scaling factors for each spectrum (all ones for the OSC)
- lightcurve_fitting.speccal.readfitsspec(filename, header=False, ext=None)[source]¶
Read a spectrum from a FITS file
- Parameters:
- filenamestr
Filename from which to read the spectrum
- headerbool, optional
If True, also return the FITS header
- extint, str, optional
FITS extension number or name in which the spectrum is stored
- Returns:
- wlarray-like
Wavelengths, typically in ångströms
- fluxarray-like
Observed fluxes in erg / (s cm2 angstrom), if units are identifiable
- hdrastropy.io.fits.header.Header, optional
FITS header, returned if
header=True
- lightcurve_fitting.speccal.readspec(f, verbose=False, return_header=False)[source]¶
Read a spectrum from a FITS or ASCII file and try to identify where and when it was observed
- Parameters:
- fstr
Filename from which to read the spectrum
- verbosebool, optional
If True, print the date and filename of the spectrum
- return_headerbool, optional
If True, return the full header as an additional return value
- Returns:
- xarray-like
Wavelengths, typically in ångströms
- yarray-like
Fluxes in erg / (s cm2 angstrom), if units are identifiable
- dateastropy.time.Time
Time at which the spectrum was observed, if identifiable (otherwise
None
)- telescopestr
Name of the telescope at which the spectrum was observed, if identifiable (otherwise
''
)- instrument: str
Name of the instrument used to observe the spectrum, if identifiable (otherwise
''
)- hdrdict-like
If
return_header=True
, metadata is returned as a FITS header or dictionary (for ASCII files)
- lightcurve_fitting.speccal.remove_duplicate_wcs(hdr, keep_number=0)[source]¶
Remove duplicate world coordinate system header keywords from a
astropy.io.fits.header.Header
- lightcurve_fitting.speccal.removebadcards(hdr)[source]¶
Remove problematic entries from a
astropy.io.fits.header.Header