P92

class dust_extinction.shapes.P92(BKG_amp=218.57142857142858, BKG_lambda=0.047, BKG_b=90.0, BKG_n=2.0, FUV_amp=18.545454545454547, FUV_lambda=0.07, FUV_b=4.0, FUV_n=6.5, NUV_amp=0.05961038961038961, NUV_lambda=0.22, NUV_b=-1.95, NUV_n=2.0, SIL1_amp=0.0026493506493506496, SIL1_lambda=9.7, SIL1_b=-1.95, SIL1_n=2.0, SIL2_amp=0.0026493506493506496, SIL2_lambda=18.0, SIL2_b=-1.8, SIL2_n=2.0, FIR_amp=0.015896103896103898, FIR_lambda=25.0, FIR_b=0.0, FIR_n=2.0, **kwargs)[source]

Bases: astropy.modeling.Fittable1DModel

P92 extinction model calculation

Parameters:
BKG_amp : float

background term amplitude

BKG_lambda : float

background term central wavelength

BKG_b : float

background term b coefficient

BKG_n : float

background term n coefficient [FIXED at n = 2]

FUV_amp : float

far-ultraviolet term amplitude

FUV_lambda : float

far-ultraviolet term central wavelength

FUV_b : float

far-ultraviolet term b coefficent

FUV_n : float

far-ultraviolet term n coefficient

NUV_amp : float

near-ultraviolet (2175 A) term amplitude

NUV_lambda : float

near-ultraviolet (2175 A) term central wavelength

NUV_b : float

near-ultraviolet (2175 A) term b coefficent

NUV_n : float

near-ultraviolet (2175 A) term n coefficient [FIXED at n = 2]

SIL1_amp : float

1st silicate feature (~10 micron) term amplitude

SIL1_lambda : float

1st silicate feature (~10 micron) term central wavelength

SIL1_b : float

1st silicate feature (~10 micron) term b coefficent

SIL1_n : float

1st silicate feature (~10 micron) term n coefficient [FIXED at n = 2]

SIL2_amp : float

2nd silicate feature (~18 micron) term amplitude

SIL2_lambda : float

2nd silicate feature (~18 micron) term central wavelength

SIL2_b : float

2nd silicate feature (~18 micron) term b coefficient

SIL2_n : float

2nd silicate feature (~18 micron) term n coefficient [FIXED at n = 2]

FIR_amp : float

far-infrared term amplitude

FIR_lambda : float

far-infrared term central wavelength

FIR_b : float

far-infrared term b coefficent

FIR_n : float

far-infrared term n coefficient [FIXED at n = 2]

Notes

P92 extinction model

From Pei (1992)

Applicable from the extreme UV to far-IR

Example showing a P92 curve with components identified.

import numpy as np
import matplotlib.pyplot as plt
import astropy.units as u

from dust_extinction.shapes import P92

fig, ax = plt.subplots()

# generate the curves and plot them
lam = np.logspace(-3.0, 3.0, num=1000)
x = (1.0/lam)/u.micron

ext_model = P92()
ax.plot(1/x,ext_model(x),label='total')

ext_model = P92(FUV_amp=0., NUV_amp=0.0,
                SIL1_amp=0.0, SIL2_amp=0.0, FIR_amp=0.0)
ax.plot(1./x,ext_model(x),label='BKG only')

ext_model = P92(NUV_amp=0.0,
                SIL1_amp=0.0, SIL2_amp=0.0, FIR_amp=0.0)
ax.plot(1./x,ext_model(x),label='BKG+FUV only')

ext_model = P92(FUV_amp=0.,
                SIL1_amp=0.0, SIL2_amp=0.0, FIR_amp=0.0)
ax.plot(1./x,ext_model(x),label='BKG+NUV only')

ext_model = P92(FUV_amp=0., NUV_amp=0.0,
                SIL2_amp=0.0)
ax.plot(1./x,ext_model(x),label='BKG+FIR+SIL1 only')

ext_model = P92(FUV_amp=0., NUV_amp=0.0,
                SIL1_amp=0.0)
ax.plot(1./x,ext_model(x),label='BKG+FIR+SIL2 only')

ext_model = P92(FUV_amp=0., NUV_amp=0.0,
                SIL1_amp=0.0, SIL2_amp=0.0)
ax.plot(1./x,ext_model(x),label='BKG+FIR only')

# Milky Way observed extinction as tabulated by Pei (1992)
MW_x = [0.21, 0.29, 0.45, 0.61, 0.80, 1.11, 1.43, 1.82,
        2.27, 2.50, 2.91, 3.65, 4.00, 4.17, 4.35, 4.57, 4.76,
        5.00, 5.26, 5.56, 5.88, 6.25, 6.71, 7.18, 7.60,
        8.00, 8.50, 9.00, 9.50, 10.00]
MW_x = np.array(MW_x)
MW_exvebv = [-3.02, -2.91, -2.76, -2.58, -2.23, -1.60, -0.78, 0.00,
             1.00, 1.30, 1.80, 3.10, 4.19, 4.90, 5.77, 6.57, 6.23,
             5.52, 4.90, 4.65, 4.60, 4.73, 4.99, 5.36, 5.91,
             6.55, 7.45, 8.45, 9.80, 11.30]
MW_exvebv = np.array(MW_exvebv)
Rv = 3.08
MW_axav = MW_exvebv/Rv + 1.0
ax.plot(1./MW_x, MW_axav, 'o', label='MW Observed')

ax.set_xscale('log')
ax.set_yscale('log')

ax.set_ylim(1e-3,10.)

ax.set_xlabel(r'$\lambda$ [$\mu$m]')
ax.set_ylabel(r'$A(x)/A(V)$')

ax.legend(loc='best')
plt.show()

(Source code, png, hires.png, pdf)

../_images/dust_extinction-shapes-P92-1.png

Attributes Summary

AbAv
BKG_amp BKG term: amplitude
BKG_b BKG term: b coefficient
BKG_lambda BKG term: center wavelength
BKG_n BKG term: n coefficient
FIR_amp FIR term: amplitude
FIR_b FIR term: b coefficient
FIR_lambda FIR term: center wavelength
FIR_n FIR term: n coefficient
FUV_amp FUV term: amplitude
FUV_b FUV term: b coefficient
FUV_lambda FUV term: center wavelength
FUV_n FUV term: n coefficient
NUV_amp NUV term: amplitude
NUV_b NUV term: b coefficient
NUV_lambda NUV term: center wavelength
NUV_n NUV term: n coefficient
SIL1_amp SIL1 term: amplitude
SIL1_b SIL1 term: b coefficient
SIL1_lambda SIL1 term: center wavelength
SIL1_n SIL1 term: n coefficient
SIL2_amp SIL2 term: amplitude
SIL2_b SIL2 term: b coefficient
SIL2_lambda SIL2 term: center wavelength
SIL2_n SIL2 term: n coefficient
fit_deriv
inputs
outputs
param_names
x_range

Methods Summary

__call__(x[, model_set_axis, …]) Evaluate this model using the given input(s) and the parameter values that were specified when the model was instantiated.
evaluate(in_x, BKG_amp, BKG_lambda, BKG_b, …) P92 function

Attributes Documentation

AbAv = 1.3246753246753247
BKG_amp

BKG term: amplitude

BKG_b

BKG term: b coefficient

BKG_lambda

BKG term: center wavelength

BKG_n

BKG term: n coefficient

FIR_amp

FIR term: amplitude

FIR_b

FIR term: b coefficient

FIR_lambda

FIR term: center wavelength

FIR_n

FIR term: n coefficient

FUV_amp

FUV term: amplitude

FUV_b

FUV term: b coefficient

FUV_lambda

FUV term: center wavelength

FUV_n

FUV term: n coefficient

NUV_amp

NUV term: amplitude

NUV_b

NUV term: b coefficient

NUV_lambda

NUV term: center wavelength

NUV_n

NUV term: n coefficient

SIL1_amp

SIL1 term: amplitude

SIL1_b

SIL1 term: b coefficient

SIL1_lambda

SIL1 term: center wavelength

SIL1_n

SIL1 term: n coefficient

SIL2_amp

SIL2 term: amplitude

SIL2_b

SIL2 term: b coefficient

SIL2_lambda

SIL2 term: center wavelength

SIL2_n

SIL2 term: n coefficient

fit_deriv = None
inputs = ('x',)
outputs = ('axav',)
param_names = ('BKG_amp', 'BKG_lambda', 'BKG_b', 'BKG_n', 'FUV_amp', 'FUV_lambda', 'FUV_b', 'FUV_n', 'NUV_amp', 'NUV_lambda', 'NUV_b', 'NUV_n', 'SIL1_amp', 'SIL1_lambda', 'SIL1_b', 'SIL1_n', 'SIL2_amp', 'SIL2_lambda', 'SIL2_b', 'SIL2_n', 'FIR_amp', 'FIR_lambda', 'FIR_b', 'FIR_n')
x_range = [0.001, 1000.0]

Methods Documentation

__call__(x, model_set_axis=None, with_bounding_box=False, fill_value=nan, equivalencies=None)

Evaluate this model using the given input(s) and the parameter values that were specified when the model was instantiated.

evaluate(in_x, BKG_amp, BKG_lambda, BKG_b, BKG_n, FUV_amp, FUV_lambda, FUV_b, FUV_n, NUV_amp, NUV_lambda, NUV_b, NUV_n, SIL1_amp, SIL1_lambda, SIL1_b, SIL1_n, SIL2_amp, SIL2_lambda, SIL2_b, SIL2_n, FIR_amp, FIR_lambda, FIR_b, FIR_n)[source]

P92 function

Parameters:
in_x: float

expects either x in units of wavelengths or frequency or assumes wavelengths in wavenumbers [1/micron]

internally wavenumbers are used

Returns:
axav: np array (float)

A(x)/A(V) extinction curve [mag]

Raises:
ValueError

Input x values outside of defined range