WG00

class dust_attenuation.radiative_transfer.WG00(tau_V, geometry='dusty', dust_type='mw', dust_distribution='clumpy')[source]

Bases: BaseAtttauVModel

Attenuation curve of Witt & Gordon (2000)

Parameters:
tau_v: float

optical depth in V band

Raises:
InputParameterError

Input Av values outside of defined range

Notes

From Witt & Gordon (2000, ApJ, Volume 528, pp. 799-816)

Example:

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

from dust_attenuation.radiative_transfer import WG00

fig, ax = plt.subplots(1,2, figsize=(10,6))

# generate the curves and plot them
# Use 1/microns for a better sampling
x = np.arange(0.35,10.0,0.1)/u.micron

x_Vband = 0.55 # microns

tau_Vs = [0.25,0.4,1.1,17.0,46.0]
for tau_V in tau_Vs[::-1]:
   att_model = WG00(tau_V = tau_V, geometry = 'cloudy',
                    dust_type = 'mw', dust_distribution = 'clumpy')
   ax[0].plot(x,att_model(1/x),label=r'$\tau_V$ = %.2f mag' % (tau_V))
   ax[1].plot(x,att_model(1/x)/att_model(x_Vband),
              label=r'$\tau_V$ = %.2f mag' % (tau_V))

ax[0].set_xlabel(r'$x$ [$\mu m^{-1}$]')
ax[0].set_ylabel(r'$Att(x)$ [mag]')
ax[1].set_xlabel(r'$x$ [$\mu m^{-1}$]')
ax[1].set_ylabel(r'$Att(x)/Att_V$')

ax[0].legend(loc='best')
ax[1].legend(loc='best')
fig.suptitle(r'CLOUDY / MW / clumpy model',size=15)
plt.tight_layout()
fig.subplots_adjust(top=0.88)

plt.show()

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

../_images/dust_attenuation-radiative_transfer-WG00-1.png

Load the attenuation curves for a given geometry, dust type and dust distribution.

Parameters:
tau_V: float

optical depth in V band

geometry: string

‘shell’, ‘cloudy’ or ‘dusty’

dust_type: string

‘mw’ or ‘smc’

dust_distribution: string

‘homogeneous’ or ‘clumpy’

Returns:
Attx: np array (float)

Att(x) attenuation curve [mag]

Attributes Summary

tau_V_range

x_range

Methods Summary

evaluate(x, tau_V)

WG00 function

get_albedo(x)

Return the albedo in function of wavelength for the corresponding dust type (SMC or MW).

get_extinction(x, tau_V)

Return the extinction at a given wavelength and V-band optical depth.

get_fdir(x, tau_V)

Return the direct attenuated stellar flux fraction at a given wavelength and V-band optical depth.

get_fesc(x, tau_V)

Return the total escaping flux fraction at a given wavelength and V-band optical depth.

get_fsca(x, tau_V)

Return the scattered flux fraction at a given wavelength and V-band optical depth.

get_scattering_phase_function(x)

Return the scattering phase function in function of wavelength for the corresponding dust type (SMC or MW).

Attributes Documentation

tau_V_range = [0.25, 50.0]
x_range = [0.1, 3.0001]

Methods Documentation

evaluate(x, tau_V)[source]

WG00 function

Parameters:
x: float

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

internally microns are used

tau_V: float

optical depth in V band

Returns:
Attx: np array (float)

Att(x) attenuation curve [mag]

Raises:
ValueError

Input x values outside of defined range

get_albedo(x)[source]

Return the albedo in function of wavelength for the corresponding dust type (SMC or MW). The albedo gives the probability a photon is scattered from a dust grain.

Parameters:
x: float

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

internally microns are used

Returns:
albedo: np array (float)

alb(x) albedo

Raises:
ValueError

Input x values outside of defined range

get_extinction(x, tau_V)[source]

Return the extinction at a given wavelength and V-band optical depth.

Parameters:
x: float

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

internally microns are used

tau_V: float

optical depth in V band

Returns:
ext: np array (float)

ext(x) extinction curve [mag]

Raises:
ValueError

Input x values outside of defined range

get_fdir(x, tau_V)[source]

Return the direct attenuated stellar flux fraction at a given wavelength and V-band optical depth.

Parameters:
x: float

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

internally microns are used

tau_V: float

optical depth in V band

Returns:
fsca: np array (float)

fsca(x) scattered flux fraction

Raises:
ValueError

Input x values outside of defined range

get_fesc(x, tau_V)[source]

Return the total escaping flux fraction at a given wavelength and V-band optical depth.

Parameters:
x: float

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

internally microns are used

tau_V: float

optical depth in V band

Returns:
fsca: np array (float)

fsca(x) scattered flux fraction

Raises:
ValueError

Input x values outside of defined range

get_fsca(x, tau_V)[source]

Return the scattered flux fraction at a given wavelength and V-band optical depth.

Parameters:
x: float

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

internally microns are used

tau_V: float

optical depth in V band

Returns:
fsca: np array (float)

fsca(x) scattered flux fraction

Raises:
ValueError

Input x values outside of defined range

get_scattering_phase_function(x)[source]

Return the scattering phase function in function of wavelength for the corresponding dust type (SMC or MW). The scattering phase function gives the angle at which the photon scatters.

Parameters:
x: float

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

internally microns are used

Returns:
g: np array (float)

g(x) scattering phase function

Raises:
ValueError

Input x values outside of defined range