Units Module (sbpy.units)

Introduction

units provides common planetary astronomy units and unit conversions, including Vega-based magnitudes. sbpy units may be added to the top-level astropy.units namespace via:

>>> import sbpy.units
>>> sbpy.units.enable()    

This will make them searchable via strings, such as u.Unit('mag(VEGA)').

Spectral Gradient Units

Spectral gradients are commonly expressed as % per 100 nm. This is too subtle a concept for astropy’s units:

>>> import astropy.units as u
>>> print(u.percent / (100 * u.nm))
0.01 % / nm
>>> print(u.Unit('% / (100 nm)'))    
ValueError: '% / (100 nm)' did not parse as unit: Syntax error parsing unit '% / 100 nm'

As a convenience, sbpy defines the hundred_nm unit, which has an appropriate string representation:

>>> from sbpy.units import hundred_nm
>>> print(u.percent / hundred_nm)
% / 100 nm

Vega Magnitudes

With the synphot package, sbpy has the capability to convert between flux densities and Vega-based magnitude systems. The conversions require a spectrum of Vega, which is provided by sbpy’s calibration system.

sbpy defines two new spectral flux density units: VEGA and JM. VEGA represents the flux density of Vega. JM represents the flux density zeropoint of the Johnson-Morgan system, assuming Vega has a magnitude of 0.03 at all wavelengths (Johnson et al. 1966, Bessell & Murphy 2012). Two magnitude units are also defined: VEGAmag and JMmag.

Unit conversions between flux density and Vega-based magnitudes use the astropy.units equivalency system. sbpy’s spectral_density_vega() provides the equivalencies, which astropy would use to convert the units. The function requires a wavelength, frequency, or bandpass:

>>> import astropy.units as u
>>> from sbpy.units import VEGAmag, spectral_density_vega
>>> wave = 5500 * u.AA
>>> m = 0 * VEGAmag
>>> fluxd = m.to('erg/(cm2 s AA)', spectral_density_vega(wave))
>>> fluxd.value                                  
3.5469235179497687e-09
>>> m = fluxd.to(VEGAmag, spectral_density_vega(wave))
>>> m.value                                      
0.0

To use a bandpass, define and pass a synphot.spectrum.SpectralElement. A limited set of bandpasses are distributed with sbpy (see Filter Bandpasses):

>>> from sbpy.units import VEGAmag, spectral_density_vega
>>> from sbpy.photometry import bandpass
>>> V = bandpass('Johnson V')
>>> m = 0.0 * VEGAmag
>>> fluxd = m.to('erg/(cm2 s AA)', spectral_density_vega(V))
>>> fluxd.value                   
3.588524658721229e-09

Reflectance Units

sbpy.units supports the conversions between reflected flux (often expressed in magnitude), and bidirectional reflectance (in unit of 1/sr) and scattering cross-section through function reflectance.

For example, the absolute magnitude of Ceres in V-band is 3.4 in VEGAmag system, the radius is 460 km, its disk-averaged bidirectional reflectance at 0 phase angle can be calculated:

>>> import numpy as np
>>> from astropy import units as u
>>> from sbpy.calib import solar_fluxd, vega_fluxd
>>> from sbpy.units import reflectance, VEGAmag, spectral_density_vega
>>>
>>> solar_fluxd.set({'V': -26.775 * VEGAmag})
<ScienceState solar_fluxd: {'V': <Magnitude -26.775 mag(VEGA)>}>
>>> vega_fluxd.set({'V': 3.5885e-08 * u.Unit('W / (m2 um)')})
<ScienceState vega_fluxd: {'V': <Quantity 3.5885e-08 W / (m2 um)>}>
>>> mag = 3.4 * VEGAmag
>>> radius = 460 * u.km
>>> cross_sec = np.pi * (radius)**2
>>> ref = mag.to('1/sr', reflectance('V', cross_section=cross_sec))
>>> print('{0:.4f}'.format(ref))
0.0287 1 / sr

reflectance works with sbpy’s spectral calibration system (see Spectral Standards and Photometric Calibration (sbpy.calib)):

>>> from sbpy.photometry import bandpass
>>> V = bandpass('Johnson V')
>>> ref = 0.0287 / u.sr
>>> cross_sec = mag.to('km2', reflectance(V, reflectance=ref))
>>> radius = np.sqrt(cross_sec / np.pi)
>>> print('{0:.2f}'.format(radius))
459.69 km

reflectance also supports conversion between a flux spectrum and a reflectance spectrum:

>>> wave = [1.046, 1.179, 1.384, 1.739, 2.416] * u.um
>>> flux = [1.636e-18, 1.157e-18, 8.523e-19, 5.262e-19, 1.9645e-19] \
...     * u.Unit('W/(m2 um)')
>>> xsec = 0.0251 * u.km**2
>>> ref = flux.to('1/sr', reflectance(wave, cross_section=xsec))
>>> print(ref)  
[0.0021763  0.00201223 0.0022041  0.00269637 0.00292785] 1 / sr

reflectance also supports dimensionless logarithmic unit for solar flux, which can be specified with set:

>>> with solar_fluxd.set({'V': -26.77 * u.mag}):
...     mag = 3.4 * u.mag
...     ref = mag.to('1/sr', reflectance('V', cross_section =
...         np.pi*(460*u.km)**2))
>>> print(ref)  
0.028786262941247264 1 / sr

Projected Sizes

With the projected_size equivalencies, one may convert between angles and lengths, for a given distance:

>>> import astropy.units as u
>>> import sbpy.units as sbu
>>> (1 * u.arcsec).to('km', sbu.projected_size(1 * u.au))
... 
<Quantity [725.27094381] km>
>>> (725.27 * u.km).to('arcsec', sbu.projected_size(1 * u.au))
... 
<Quantity [1.00] arcsec>

Reference/API

Common planetary astronomy units.

Available Units

Unit

Description

Represents

Aliases

SI Prefixes

JM

Johnson-Morgan magnitude system flux density zeropoint (Johnson et al 1966; Bessell & Murphy 2012).

\(\mathrm{1.0280163\,VEGA}\)

JMflux

No

VEGA

Spectral flux density of Vega.

VEGAflux

No

Functions

enable()

Enable sbpy units in the top-level astropy.units namespace.

projected_size(eph)

Angular size and projected linear size equivalency.

reflectance(wfb[, cross_section, reflectance])

Reflectance related equivalencies.

spectral_density_vega(wfb)

Flux density equivalencies with Vega-based magnitude systems.

Variables

JM

Johnson-Morgan magnitude system flux density zeropoint (Johnson et al 1966; Bessell & Murphy 2012).

JMmag

Johnson-Morgan magnitude system: Vega is 0.03 mag at all wavelengths (Johnson et al 1966; Bessell & Murphy 2012).

VEGA

Spectral flux density of Vega.

VEGAmag

Vega-based magnitude: Vega is 0 mag at all wavelengths

hundred_nm

Convenience unit for expressing spectral gradients.