Interstellar Dust Extinction¶
dust_extinction
is a python package to provide models of interstellar dust
extinction curves.
Extinction describes the effects of dust on observations of single star due to the dust along the line-of-sight to a star removing flux by absorbing photons and scattering photons out of the line-of-sight. The wavelength dependence of dust extinction (also know as extinction curves) provides fundamental information about the size, composition, and shape of interstellar dust grain. In general, extinction models are used to model or correct the effects of dust on observations.
In contrast, dust attenuation refers to the effects of dust on the measurements of groups of stars mixed with dust or stars with circumstellar dust. See Extinction versus Attenuation. For attenuation models, see the dust_attenuation package.
This package is an astropy affiliated package and uses the astropy.modeling framework.
User Documentation¶
Dev Documentation¶
Installation¶
Repository¶
GitHub: dust_extinction
Quick Start¶
Extinction Curve¶
How to get the A(x)/A(V) extinction curve for an input set of x wavelength values.
Define a model, specifically the G23 model with an R(V) = 3.1.
from dust_extinction.parameter_averages import G23
# define the model
extmod = G23(Rv=3.1)
Define the wavelengths
wavelengths = np.logspace(np.log10(0.1), np.log10(30.0), num=1000)*u.micron
Get the extinction values in A(lambda)/A(V) units.
ext = extmod(wavelengths)
(Source code
, png
, hires.png
, pdf
)
Extinguish/Unextinguish a Spectrum¶
How to extinguish (redden) and unextinguish (deredden) a spectrum:
Generate a spectrum to use. In this case a blackbody model, but can be an
observed spectrum. The dust_extinction
models are unit aware and the
wavelength array should have astropy.units associated with it.
import numpy as np
from astropy.modeling.models import BlackBody
import astropy.units as u
# wavelengths and spectrum are 1D arrays
# wavelengths between 1000 and 30000 A
wavelengths = np.logspace(np.log10(0.1), np.log10(30.0), num=1000)*u.micron
bb_lam = BlackBody(10000*u.K, scale=1.0 * u.erg / (u.cm ** 2 * u.AA * u.s * u.sr))
spectrum = bb_lam(wavelengths)
Define a model, specifically the G23 model with an R(V) = 3.1.
from dust_extinction.parameter_averages import G23
# define the model
ext = G23(Rv=3.1)
Extinguish (redden) a spectrum with a screen of G23 dust with an E(B-V) of 0.5. Can also specify the dust column with Av (this case equivalent to Av = 0.5*Rv = 1.55).
# extinguish (redden) the spectrum
spectrum_ext = spectrum*ext.extinguish(wavelengths, Ebv=0.5)
Unextinguish (deredden) a spectrum with a screen of G23 dust with the equivalent A(V) column.
# unextinguish (deredden) the spectrum
# Av = 1.55 = R(V) * E(B-V) = 3.1 * 0.5
spectrum_noext = spectrum_ext/ext.extinguish(wavelengths, Av=1.55)
(Source code
, png
, hires.png
, pdf
)
Reporting Issues¶
If you have found a bug in dust_extinction
please report it by creating a
new issue on the dust_extinction
GitHub issue tracker.
Please include an example that demonstrates the issue sufficiently so that the developers can reproduce and fix the problem. You may also be asked to provide information about your operating system and a full Python stack trace. The developers will walk you through obtaining a stack trace if it is necessary.
Contributing¶
Like the Astropy project, dust_extinction
is made both by and for its
users. We accept contributions at all levels, spanning the gamut from fixing a
typo in the documentation to developing a major new feature. We welcome
contributors who will abide by the Python Software Foundation Code of Conduct.
dust_extinction
follows the same workflow and coding guidelines as
Astropy. Take a look at the astropy
developer documentation for
guidelines.
For the complete list of contributors please see the dust_extinction contributors page on Github.