trident.absorption_spectrum.absorption_spectrum.AbsorptionSpectrum
- class trident.absorption_spectrum.absorption_spectrum.AbsorptionSpectrum(lambda_min, lambda_max, n_lambda=None, dlambda=None, bin_space='wavelength')[source]
Base class for generating absorption spectra. This code was originally based in yt and more restrictive in terms of what development was allowed, so the
SpectrumGenerator
subclass has more advanced functionality built on top of this. The base algorithm and functionality for spectral generation occurs here though.Note
The preferred method for generating spectra is using
SpectrumGenerator
.Parameters
- Lambda_min
float, YTQuantity, or ‘auto’
lower wavelength bound in angstroms or velocity bound in km/s (if bin_space set to ‘velocity’). If set to ‘auto’, the lower bound will be automatically adjusted to encompass all absorption lines. The window will not be expanded for continuum features, only absorption lines.
- Lambda_max
float, YTQuantity, or ‘auto’
upper wavelength bound in angstroms or velocity bound in km/s (if bin_space set to ‘velocity’). If set to ‘auto’, the upper bound will be automatically adjusted to encompass all absorption lines. The window will not be expanded for continuum features, only absorption lines.
- N_lambda
optional, int
number of bins. This cannot be set when setting either lambda_min or lambda_max to auto.
- Dlambda
optional, float or YTQuantity
size of the wavelength bins in angstroms or velocity bins in km/s.
- Bin_space
‘wavelength’ or ‘velocity’
Sets the dimension in which spectra are created. If set to wavelength, the resulting spectra are flux (or tau) vs. observed wavelength. If set to velocity, the spectra are flux vs. velocity offset from the rest wavelength of the absorption line. Default: wavelength
Methods
Add a continuum feature that follows a power-law. |
|
Add an absorption line to the list of lines included in the spectrum. |
|
Approximate the flux error for a spectrum. |
|
Make spectrum from ray data using the line list. |
Attributes
This is the optical depth array for the current absorption line being deposited. |
|
The lambda field. |
|
This is the total optical depth of all lines and continua. |