Units and Quantities (astropy.units)

Introduction

astropy.units handles defining, converting between, and performing arithmetic with physical quantities, such as meters, seconds, Hz, etc. It also handles logarithmic units such as magnitude and decibel.

astropy.units does not know spherical geometry or sexagesimal (hours, min, sec): if you want to deal with celestial coordinates, see the astropy.coordinates package.

Getting Started

Most users of the astropy.units package will work with Quantity objects: the combination of a value and a unit. The most convenient way to create a Quantity is to multiply or divide a value by one of the built-in units. It works with scalars, sequences, and numpy arrays.

Examples

To create a Quantity object:

>>> from astropy import units as u
>>> 42.0 * u.meter  
<Quantity  42. m>
>>> [1., 2., 3.] * u.m  
<Quantity [1., 2., 3.] m>
>>> import numpy as np
>>> np.array([1., 2., 3.]) * u.m  
<Quantity [1., 2., 3.] m>

You can get the unit and value from a Quantity using the unit and value members:

>>> q = 42.0 * u.meter
>>> q.value
42.0
>>> q.unit
Unit("m")

From this basic building block, it is possible to start combining quantities with different units:

>>> 15.1 * u.meter / (32.0 * u.second)  
<Quantity 0.471875 m / s>
>>> 3.0 * u.kilometer / (130.51 * u.meter / u.second)  
<Quantity 0.022986744310780783 km s / m>
>>> (3.0 * u.kilometer / (130.51 * u.meter / u.second)).decompose()  
<Quantity 22.986744310780782 s>

Unit conversion is done using the to() method, which returns a new Quantity in the given unit:

>>> x = 1.0 * u.parsec
>>> x.to(u.km)  
<Quantity 30856775814671.914 km>

It is also possible to work directly with units at a lower level, for example, to create custom units:

>>> from astropy.units import imperial

>>> cms = u.cm / u.s
>>> # ...and then use some imperial units
>>> mph = imperial.mile / u.hour

>>> # And do some conversions
>>> q = 42.0 * cms
>>> q.to(mph)  
<Quantity 0.939513242662849 mi / h>

Units that “cancel out” become a special unit called the “dimensionless unit”:

>>> u.m / u.m
Unit(dimensionless)

To create a basic dimensionless quantity, multiply a value by the unscaled dimensionless unit:

>>> q = 1.0 * u.dimensionless_unscaled
>>> q.unit
Unit(dimensionless)

astropy.units is able to match compound units against the units it already knows about:

>>> (u.s ** -1).compose()  
[Unit("Bq"), Unit("Hz"), Unit("2.7027e-11 Ci")]

And it can convert between unit systems, such as SI or CGS:

>>> (1.0 * u.Pa).cgs
<Quantity 10. P / s>

The units mag, dex, and dB are special, being logarithmic units, for which a value is the logarithm of a physical quantity in a given unit. These can be used with a physical unit in parentheses to create a corresponding logarithmic quantity:

>>> -2.5 * u.mag(u.ct / u.s)
<Magnitude -2.5 mag(ct / s)>
>>> from astropy import constants as c
>>> u.Dex((c.G * u.M_sun / u.R_sun**2).cgs)  
<Dex 4.438067627303133 dex(cm / s2)>

astropy.units also handles equivalencies, such as that between wavelength and frequency. To use that feature, equivalence objects are passed to the to() conversion method. For instance, a conversion from wavelength to frequency does not normally work:

>>> (1000 * u.nm).to(u.Hz)  
Traceback (most recent call last):
  ...
UnitConversionError: 'nm' (length) and 'Hz' (frequency) are not convertible

But by passing an equivalency list, in this case spectral(), it does:

>>> (1000 * u.nm).to(u.Hz, equivalencies=u.spectral())  
<Quantity  2.99792458e+14 Hz>

Quantities and units can be printed nicely to strings using the Format String Syntax. Format specifiers (like 0.03f) in strings will be used to format the quantity value:

>>> q = 15.1 * u.meter / (32.0 * u.second)
>>> q  
<Quantity 0.471875 m / s>
>>> f"{q:0.03f}"
'0.472 m / s'

The value and unit can also be formatted separately. Format specifiers for units can be used to choose the unit formatter:

>>> q = 15.1 * u.meter / (32.0 * u.second)
>>> q  
<Quantity 0.471875 m / s>
>>> f"{q.value:0.03f} {q.unit:FITS}"
'0.472 m s-1'

Using astropy.units

Acknowledgments

This code was originally based on the pynbody units module written by Andrew Pontzen, who has granted the Astropy Project permission to use the code under a BSD license.

See Also

Performance Tips

If you are attaching units to arrays to make Quantity objects, multiplying arrays by units will result in the array being copied in memory, which will slow things down. Furthermore, if you are multiplying an array by a composite unit, the array will be copied for each individual multiplication. Thus, in the following case, the array is copied four successive times:

In [1]: import numpy as np

In [2]: from astropy import units as u

In [3]: rng = np.random.default_rng()

In [4]: array = rng.random(10000000)

In [5]: %timeit array * u.m / u.s / u.kg / u.sr
92.5 ms ± 2.52 ms per loop (mean ± std. dev. of 7 runs, 10 loops each)

There are several ways to speed this up. First, when you are using composite units, ensure that the entire unit gets evaluated first, then attached to the array. You can do this by using parentheses as for any other operation:

In [6]: %timeit array * (u.m / u.s / u.kg / u.sr)
21.5 ms ± 886 µs per loop (mean ± std. dev. of 7 runs, 10 loops each)

In this case, this has sped things up by a factor of 4. If you use a composite unit several times in your code then you can define a variable for it:

In [7]: UNIT_MSKGSR = u.m / u.s / u.kg / u.sr

In [8]: %timeit array * UNIT_MSKGSR
22.2 ms ± 551 µs per loop (mean ± std. dev. of 7 runs, 10 loops each)

In this case and the case with brackets, the array is still copied once when creating the Quantity. If you want to avoid any copies altogether, you can make use of the << operator to attach the unit to the array:

In [9]: %timeit array << u.m / u.s / u.kg / u.sr
47.1 µs ± 5.77 µs per loop (mean ± std. dev. of 7 runs, 10000 loops each)

Note that these are now microseconds, so this is 2000x faster than the original case with no brackets. Note that brackets are not needed when using << since * and / have a higher precedence, so the unit will be evaluated first. When using <<, be aware that because the data is not being copied, changing the original array will also change the Quantity object.

Note that for composite units, you will definitely see an impact if you can pre-compute the composite unit:

In [10]: %timeit array << UNIT_MSKGSR
6.51 µs ± 112 ns per loop (mean ± std. dev. of 7 runs, 100000 loops each)

Which is over 10000x faster than the original example. See Creating and Converting Quantities without Copies for more details about the << operator.

Reference/API

astropy.units.quantity Module

This module defines the Quantity object, which represents a number with some associated units. Quantity objects support operations like ordinary numbers, but will deal with unit conversions internally.

Functions

allclose(a, b[, rtol, atol, equal_nan])

Whether two arrays are element-wise equal within a tolerance.

isclose(a, b[, rtol, atol, equal_nan])

Return a boolean array where two arrays are element-wise equal within a tolerance.

Classes

Quantity(value[, unit, dtype, copy, order, ...])

A Quantity represents a number with some associated unit.

SpecificTypeQuantity(value[, unit, dtype, ...])

Superclass for Quantities of specific physical type.

QuantityInfoBase([bound])

QuantityInfo([bound])

Container for meta information like name, description, format.

Class Inheritance Diagram

Inheritance diagram of astropy.units.quantity.Quantity, astropy.units.quantity.SpecificTypeQuantity, astropy.units.quantity.QuantityInfoBase, astropy.units.quantity.QuantityInfo

astropy.units Package

This subpackage contains classes and functions for defining and converting between different physical units.

This code is adapted from the pynbody units module written by Andrew Pontzen, who has granted the Astropy project permission to use the code under a BSD license.

Functions

add_enabled_aliases(aliases)

Add aliases for units.

add_enabled_equivalencies(equivalencies)

Adds to the equivalencies enabled in the unit registry.

add_enabled_units(units)

Adds to the set of units enabled in the unit registry.

allclose(a, b[, rtol, atol, equal_nan])

Whether two arrays are element-wise equal within a tolerance.

beam_angular_area(beam_area)

Convert between the beam unit, which is commonly used to express the area of a radio telescope resolution element, and an area on the sky.

brightness_temperature(frequency[, beam_area])

Defines the conversion between Jy/sr and "brightness temperature", \(T_B\), in Kelvins.

def_physical_type(unit, name)

Add a mapping between a unit and the corresponding physical type(s).

def_unit(s[, represents, doc, format, ...])

Factory function for defining new units.

dimensionless_angles()

Allow angles to be equivalent to dimensionless (with 1 rad = 1 m/m = 1).

doppler_optical(rest)

Return the equivalency pairs for the optical convention for velocity.

doppler_radio(rest)

Return the equivalency pairs for the radio convention for velocity.

doppler_redshift()

Returns the equivalence between Doppler redshift (unitless) and radial velocity.

doppler_relativistic(rest)

Return the equivalency pairs for the relativistic convention for velocity.

get_current_unit_registry()

get_physical_type(obj)

Return the physical type that corresponds to a unit (or another physical type representation).

isclose(a, b[, rtol, atol, equal_nan])

Return a boolean array where two arrays are element-wise equal within a tolerance.

logarithmic()

Allow logarithmic units to be converted to dimensionless fractions

mass_energy()

Returns a list of equivalence pairs that handle the conversion between mass and energy.

molar_mass_amu()

Returns the equivalence between amu and molar mass.

parallax()

Returns a list of equivalence pairs that handle the conversion between parallax angle and distance.

pixel_scale(pixscale)

Convert between pixel distances (in units of pix) and other units, given a particular pixscale.

plate_scale(platescale)

Convert between lengths (to be interpreted as lengths in the focal plane) and angular units with a specified platescale.

quantity_input([func])

A decorator for validating the units of arguments to functions.

set_enabled_aliases(aliases)

Set aliases for units.

set_enabled_equivalencies(equivalencies)

Sets the equivalencies enabled in the unit registry.

set_enabled_units(units)

Sets the units enabled in the unit registry.

spectral()

Returns a list of equivalence pairs that handle spectral wavelength, wave number, frequency, and energy equivalencies.

spectral_density(wav[, factor])

Returns a list of equivalence pairs that handle spectral density with regard to wavelength and frequency.

temperature()

Convert between Kelvin, Celsius, Rankine and Fahrenheit here because Unit and CompositeUnit cannot do addition or subtraction properly.

temperature_energy()

Convert between Kelvin and keV(eV) to an equivalent amount.

thermodynamic_temperature(frequency[, T_cmb])

Defines the conversion between Jy/sr and "thermodynamic temperature", \(T_{CMB}\), in Kelvins.

zero_point_flux(flux0)

An equivalency for converting linear flux units ("maggys") defined relative to a standard source into a standardized system.

Classes

CompositeUnit(scale, bases, powers[, ...])

Create a composite unit using expressions of previously defined units.

Decibel(value[, unit, dtype, copy, order, ...])

DecibelUnit([physical_unit, function_unit])

Logarithmic physical units expressed in dB

Dex(value[, unit, dtype, copy, order, ...])

DexUnit([physical_unit, function_unit])

Logarithmic physical units expressed in magnitudes

Equivalency(equiv_list[, name, kwargs])

A container for a units equivalency.

FunctionQuantity(value[, unit, dtype, copy, ...])

A representation of a (scaled) function of a number with a unit.

FunctionUnitBase([physical_unit, function_unit])

Abstract base class for function units.

IrreducibleUnit(st[, doc, format, namespace])

Irreducible units are the units that all other units are defined in terms of.

LogQuantity(value[, unit, dtype, copy, ...])

A representation of a (scaled) logarithm of a number with a unit

LogUnit([physical_unit, function_unit])

Logarithmic unit containing a physical one

MagUnit([physical_unit, function_unit])

Logarithmic physical units expressed in magnitudes

Magnitude(value[, unit, dtype, copy, order, ...])

NamedUnit(st[, doc, format, namespace])

The base class of units that have a name.

PhysicalType(unit, physical_types)

Represents the physical type(s) that are dimensionally compatible with a set of units.

PrefixUnit([s, represents, format, ...])

A unit that is simply a SI-prefixed version of another unit.

Quantity(value[, unit, dtype, copy, order, ...])

A Quantity represents a number with some associated unit.

QuantityInfo([bound])

Container for meta information like name, description, format.

QuantityInfoBase([bound])

SpecificTypeQuantity(value[, unit, dtype, ...])

Superclass for Quantities of specific physical type.

StructuredUnit(units[, names])

Container for units for a structured Quantity.

Unit([s, represents, format, namespace, ...])

The main unit class.

UnitBase()

Abstract base class for units.

UnitConversionError

Used specifically for errors related to converting between units or interpreting units in terms of other units.

UnitTypeError

Used specifically for errors in setting to units not allowed by a class.

UnitsError

The base class for unit-specific exceptions.

UnitsWarning

The base class for unit-specific warnings.

UnrecognizedUnit(st[, doc, format, namespace])

A unit that did not parse correctly.

Class Inheritance Diagram

Inheritance diagram of astropy.units.core.CompositeUnit, astropy.units.function.logarithmic.Decibel, astropy.units.function.logarithmic.DecibelUnit, astropy.units.function.logarithmic.Dex, astropy.units.function.logarithmic.DexUnit, astropy.units.equivalencies.Equivalency, astropy.units.function.core.FunctionQuantity, astropy.units.function.core.FunctionUnitBase, astropy.units.core.IrreducibleUnit, astropy.units.function.logarithmic.LogQuantity, astropy.units.function.logarithmic.LogUnit, astropy.units.function.logarithmic.MagUnit, astropy.units.function.logarithmic.Magnitude, astropy.units.core.NamedUnit, astropy.units.physical.PhysicalType, astropy.units.core.PrefixUnit, astropy.units.quantity.Quantity, astropy.units.quantity.QuantityInfo, astropy.units.quantity.QuantityInfoBase, astropy.units.quantity.SpecificTypeQuantity, astropy.units.structured.StructuredUnit, astropy.units.core.Unit, astropy.units.core.UnitBase, astropy.units.core.UnitConversionError, astropy.units.core.UnitTypeError, astropy.units.core.UnitsError, astropy.units.core.UnitsWarning, astropy.units.core.UnrecognizedUnit

astropy.units.format Package

A collection of different unit formats.

Functions

get_format([format])

Get a formatter by name.

Classes

Base(*args, **kwargs)

The abstract base class of all unit formats.

Generic(*args, **kwargs)

A "generic" format.

CDS(*args, **kwargs)

Support the Centre de Données astronomiques de Strasbourg Standards for Astronomical Catalogues 2.0 format, and the complete set of supported units.

Console(*args, **kwargs)

Output-only format for to display pretty formatting at the console.

Fits(*args, **kwargs)

The FITS standard unit format.

Latex(*args, **kwargs)

Output LaTeX to display the unit based on IAU style guidelines.

LatexInline(*args, **kwargs)

Output LaTeX to display the unit based on IAU style guidelines with negative powers.

OGIP(*args, **kwargs)

Support the units in Office of Guest Investigator Programs (OGIP) FITS files.

Unicode(*args, **kwargs)

Output-only format to display pretty formatting at the console using Unicode characters.

Unscaled(*args, **kwargs)

A format that doesn't display the scale part of the unit, other than that, it is identical to the Generic format.

VOUnit(*args, **kwargs)

The IVOA standard for units used by the VO.

Class Inheritance Diagram

Inheritance diagram of astropy.units.format.base.Base, astropy.units.format.generic.Generic, astropy.units.format.cds.CDS, astropy.units.format.console.Console, astropy.units.format.fits.Fits, astropy.units.format.latex.Latex, astropy.units.format.latex.LatexInline, astropy.units.format.ogip.OGIP, astropy.units.format.unicode_format.Unicode, astropy.units.format.generic.Unscaled, astropy.units.format.vounit.VOUnit

astropy.units.si Module

This package defines the SI units. They are also available in the astropy.units namespace.

Available Units

Unit

Description

Represents

Aliases

SI Prefixes

a

annum (a)

\(\mathrm{365.25\,d}\)

annum

Yes

A

ampere: base unit of electric current in SI

ampere, amp

Yes

Angstrom

ångström: 10 ** -10 m

\(\mathrm{0.1\,nm}\)

AA, angstrom

Yes

arcmin

arc minute: angular measurement

\(\mathrm{0.016666667\,{}^{\circ}}\)

arcminute

Yes

arcsec

arc second: angular measurement

\(\mathrm{0.00027777778\,{}^{\circ}}\)

arcsecond

Yes

Bq

becquerel: unit of radioactivity

\(\mathrm{\frac{1}{s}}\)

becquerel

No

C

coulomb: electric charge

\(\mathrm{A\,s}\)

coulomb

Yes

cd

candela: base unit of luminous intensity in SI

candela

Yes

Ci

curie: unit of radioactivity

\(\mathrm{3.7 \times 10^{10}\,Bq}\)

curie

No

d

day (d)

\(\mathrm{24\,h}\)

day

Yes

deg

degree: angular measurement 1/360 of full rotation

\(\mathrm{0.017453293\,rad}\)

degree

Yes

deg_C

Degrees Celsius

Celsius

No

eV

Electron Volt

\(\mathrm{1.6021766 \times 10^{-19}\,J}\)

electronvolt

Yes

F

Farad: electrical capacitance

\(\mathrm{\frac{C}{V}}\)

Farad, farad

Yes

fortnight

fortnight

\(\mathrm{2\,wk}\)

No

g

gram (g)

\(\mathrm{0.001\,kg}\)

gram

Yes

h

hour (h)

\(\mathrm{3600\,s}\)

hour, hr

Yes

H

Henry: inductance

\(\mathrm{\frac{Wb}{A}}\)

Henry, henry

Yes

hourangle

hour angle: angular measurement with 24 in a full circle

\(\mathrm{15\,{}^{\circ}}\)

No

Hz

Frequency

\(\mathrm{\frac{1}{s}}\)

Hertz, hertz

Yes

J

Joule: energy

\(\mathrm{N\,m}\)

Joule, joule

Yes

K

Kelvin: temperature with a null point at absolute zero.

Kelvin

Yes

kg

kilogram: base unit of mass in SI.

kilogram

No

l

liter: metric unit of volume

\(\mathrm{1000\,cm^{3}}\)

L, liter

Yes

lm

lumen: luminous flux

\(\mathrm{cd\,sr}\)

lumen

Yes

lx

lux: luminous emittance

\(\mathrm{\frac{lm}{m^{2}}}\)

lux

Yes

m

meter: base unit of length in SI

meter

Yes

mas

milli arc second: angular measurement

\(\mathrm{0.001\,{}^{\prime\prime}}\)

No

micron

micron: alias for micrometer (um)

\(\mathrm{\mu m}\)

No

min

minute (min)

\(\mathrm{60\,s}\)

minute

Yes

mol

mole: amount of a chemical substance in SI.

mole

Yes

N

Newton: force

\(\mathrm{\frac{kg\,m}{s^{2}}}\)

Newton, newton

Yes

Ohm

Ohm: electrical resistance

\(\mathrm{\frac{V}{A}}\)

ohm

Yes

Pa

Pascal: pressure

\(\mathrm{\frac{J}{m^{3}}}\)

Pascal, pascal

Yes

%

percent: one hundredth of unity, factor 0.01

\(\mathrm{0.01\,}\)

pct

No

rad

radian: angular measurement of the ratio between the length on an arc and its radius

radian

Yes

s

second: base unit of time in SI.

second

Yes

S

Siemens: electrical conductance

\(\mathrm{\frac{A}{V}}\)

Siemens, siemens

Yes

sday

Sidereal day (sday) is the time of one rotation of the Earth.

\(\mathrm{86164.091\,s}\)

No

sr

steradian: unit of solid angle in SI

\(\mathrm{rad^{2}}\)

steradian

Yes

t

Metric tonne

\(\mathrm{1000\,kg}\)

tonne

No

T

Tesla: magnetic flux density

\(\mathrm{\frac{Wb}{m^{2}}}\)

Tesla, tesla

Yes

uas

micro arc second: angular measurement

\(\mathrm{1 \times 10^{-6}\,{}^{\prime\prime}}\)

No

V

Volt: electric potential or electromotive force

\(\mathrm{\frac{J}{C}}\)

Volt, volt

Yes

W

Watt: power

\(\mathrm{\frac{J}{s}}\)

Watt, watt

Yes

Wb

Weber: magnetic flux

\(\mathrm{V\,s}\)

Weber, weber

Yes

wk

week (wk)

\(\mathrm{7\,d}\)

week

No

yr

year (yr)

\(\mathrm{365.25\,d}\)

year

Yes

astropy.units.cgs Module

This package defines the CGS units. They are also available in the top-level astropy.units namespace.

Available Units

Unit

Description

Represents

Aliases

SI Prefixes

abC

abcoulomb: CGS (EMU) of charge

\(\mathrm{Bi\,s}\)

abcoulomb

No

Ba

Barye: CGS unit of pressure

\(\mathrm{\frac{g}{cm\,s^{2}}}\)

Barye, barye

Yes

Bi

Biot: CGS (EMU) unit of current

\(\mathrm{\frac{cm^{1/2}\,g^{1/2}}{s}}\)

Biot, abA, abampere

No

C

coulomb: electric charge

\(\mathrm{A\,s}\)

coulomb

No

cd

candela: base unit of luminous intensity in SI

candela

No

cm

centimeter (cm)

\(\mathrm{cm}\)

centimeter

No

D

Debye: CGS unit of electric dipole moment

\(\mathrm{3.3333333 \times 10^{-30}\,C\,m}\)

Debye, debye

Yes

deg_C

Degrees Celsius

Celsius

No

dyn

dyne: CGS unit of force

\(\mathrm{\frac{cm\,g}{s^{2}}}\)

dyne

Yes

erg

erg: CGS unit of energy

\(\mathrm{\frac{cm^{2}\,g}{s^{2}}}\)

Yes

Fr

Franklin: CGS (ESU) unit of charge

\(\mathrm{\frac{cm^{3/2}\,g^{1/2}}{s}}\)

Franklin, statcoulomb, statC, esu

No

g

gram (g)

\(\mathrm{0.001\,kg}\)

gram

No

G

Gauss: CGS unit for magnetic field

\(\mathrm{0.0001\,T}\)

Gauss, gauss

Yes

Gal

Gal: CGS unit of acceleration

\(\mathrm{\frac{cm}{s^{2}}}\)

gal

Yes

K

Kelvin: temperature with a null point at absolute zero.

Kelvin

No

k

kayser: CGS unit of wavenumber

\(\mathrm{\frac{1}{cm}}\)

Kayser, kayser

Yes

mol

mole: amount of a chemical substance in SI.

mole

No

Mx

Maxwell: CGS unit for magnetic flux

\(\mathrm{1 \times 10^{-8}\,Wb}\)

Maxwell, maxwell

No

P

poise: CGS unit of dynamic viscosity

\(\mathrm{\frac{g}{cm\,s}}\)

poise

Yes

rad

radian: angular measurement of the ratio between the length on an arc and its radius

radian

No

s

second: base unit of time in SI.

second

No

sr

steradian: unit of solid angle in SI

\(\mathrm{rad^{2}}\)

steradian

No

St

stokes: CGS unit of kinematic viscosity

\(\mathrm{\frac{cm^{2}}{s}}\)

stokes

Yes

statA

statampere: CGS (ESU) unit of current

\(\mathrm{\frac{Fr}{s}}\)

statampere

No

astropy.units.astrophys Module

This package defines the astrophysics-specific units. They are also available in the astropy.units namespace.

Available Units

Unit

Description

Represents

Aliases

SI Prefixes

adu

adu

Yes

AU

astronomical unit: approximately the mean Earth–Sun distance.

\(\mathrm{1.4959787 \times 10^{11}\,m}\)

au, astronomical_unit

Yes

beam

beam

Yes

bin

bin

Yes

chan

chan

Yes

ct

count (ct)

count

Yes

DN

dn (DN)

dn

No

earthMass

Earth mass

\(\mathrm{5.9721679 \times 10^{24}\,kg}\)

M_earth, Mearth

No

earthRad

Earth radius

\(\mathrm{6378100\,m}\)

R_earth, Rearth

No

electron

Number of electrons

No

jupiterMass

Jupiter mass

\(\mathrm{1.8981246 \times 10^{27}\,kg}\)

M_jup, Mjup, M_jupiter, Mjupiter

No

jupiterRad

Jupiter radius

\(\mathrm{71492000\,m}\)

R_jup, Rjup, R_jupiter, Rjupiter

No

Jy

Jansky: spectral flux density

\(\mathrm{1 \times 10^{-26}\,\frac{W}{Hz\,m^{2}}}\)

Jansky, jansky

Yes

lsec

Light second

\(\mathrm{2.9979246 \times 10^{8}\,m}\)

lightsecond

No

lyr

Light year

\(\mathrm{9.4607305 \times 10^{15}\,m}\)

lightyear

Yes

pc

parsec: approximately 3.26 light-years.

\(\mathrm{3.0856776 \times 10^{16}\,m}\)

parsec

Yes

ph

photon (ph)

photon

Yes

R

Rayleigh: photon flux

\(\mathrm{7.9577472 \times 10^{8}\,\frac{ph}{s\,sr\,m^{2}}}\)

Rayleigh, rayleigh

Yes

Ry

Rydberg: Energy of a photon whose wavenumber is the Rydberg constant

\(\mathrm{13.605693\,eV}\)

rydberg

Yes

solLum

Solar luminance

\(\mathrm{3.828 \times 10^{26}\,W}\)

L_sun, Lsun

No

solMass

Solar mass

\(\mathrm{1.9884099 \times 10^{30}\,kg}\)

M_sun, Msun

No

solRad

Solar radius

\(\mathrm{6.957 \times 10^{8}\,m}\)

R_sun, Rsun

No

Sun

Sun

No

astropy.units.misc Module

This package defines miscellaneous units. They are also available in the astropy.units namespace.

Available Units

Unit

Description

Represents

Aliases

SI Prefixes

bar

bar: pressure

\(\mathrm{100000\,Pa}\)

Yes

barn

barn: unit of area used in HEP

\(\mathrm{1 \times 10^{-28}\,m^{2}}\)

barn

Yes

bit

b (bit)

b

Yes

byte

B (byte)

\(\mathrm{8\,bit}\)

B

Yes

cycle

cycle: angular measurement, a full turn or rotation

\(\mathrm{6.2831853\,rad}\)

cy

No

M_e

Electron mass

\(\mathrm{9.1093837 \times 10^{-31}\,kg}\)

No

M_p

Proton mass

\(\mathrm{1.6726219 \times 10^{-27}\,kg}\)

No

pix

pixel (pix)

pixel

Yes

spat

spat: the solid angle of the sphere, 4pi sr

\(\mathrm{12.566371\,sr}\)

sp

No

Torr

Unit of pressure based on an absolute scale, now defined as exactly 1/760 of a standard atmosphere

\(\mathrm{133.32237\,Pa}\)

torr

Yes

u

Unified atomic mass unit

\(\mathrm{1.6605391 \times 10^{-27}\,kg}\)

Da, Dalton

Yes

vox

voxel (vox)

voxel

Yes

astropy.units.function.units Module

This package defines units that can also be used as functions of other units. If called, their arguments are used to initialize the corresponding function unit (e.g., u.mag(u.ct/u.s)). Note that the prefixed versions cannot be called, as it would be unclear what, e.g., u.mmag(u.ct/u.s) would mean.

Available Units

Unit

Description

Represents

Aliases

SI Prefixes

dB

Decibel: ten per base 10 logarithmic unit

\(\mathrm{0.1\,dex}\)

decibel

No

dex

Dex: Base 10 logarithmic unit

No

mag

Astronomical magnitude: -2.5 per base 10 logarithmic unit

\(\mathrm{-0.4\,dex}\)

Yes

astropy.units.photometric Module

This module defines magnitude zero points and related photometric quantities.

The corresponding magnitudes are given in the description of each unit (the actual definitions are in logarithmic).

Available Units

Unit

Description

Represents

Aliases

SI Prefixes

AB

AB magnitude zero flux density (magnitude ABmag).

\(\mathrm{3.6307805 \times 10^{-20}\,\frac{erg}{Hz\,s\,cm^{2}}}\)

ABflux

No

Bol

Luminosity corresponding to absolute bolometric magnitude zero (magnitude M_bol).

\(\mathrm{3.0128 \times 10^{28}\,W}\)

L_bol

No

bol

Irradiance corresponding to appparent bolometric magnitude zero (magnitude m_bol).

\(\mathrm{2.3975101 \times 10^{25}\,\frac{W}{pc^{2}}}\)

f_bol

No

mgy

Maggies - a linear flux unit that is the flux for a mag=0 object.To tie this onto a specific calibrated unit system, the zero_point_flux equivalency should be used.

maggy

Yes

ST

ST magnitude zero flux density (magnitude STmag).

\(\mathrm{3.6307805 \times 10^{-9}\,\frac{erg}{\mathring{A}\,s\,cm^{2}}}\)

STflux

No

Functions

zero_point_flux(flux0)

An equivalency for converting linear flux units ("maggys") defined relative to a standard source into a standardized system.

astropy.units.imperial Module

This package defines colloquially used Imperial units. They are available in the astropy.units.imperial namespace, but not in the top-level astropy.units namespace, e.g.:

>>> import astropy.units as u
>>> mph = u.imperial.mile / u.hour
>>> mph
Unit("mi / h")

To include them in compose and the results of find_equivalent_units, do:

>>> import astropy.units as u
>>> u.imperial.enable()  
Available Units

Unit

Description

Represents

Aliases

SI Prefixes

ac

International acre

\(\mathrm{43560\,ft^{2}}\)

acre

No

BTU

British thermal unit

\(\mathrm{1.0550559\,kJ}\)

btu

No

cal

Thermochemical calorie: pre-SI metric unit of energy

\(\mathrm{4.184\,J}\)

calorie

No

cup

U.S.

\(\mathrm{0.5\,pint}\)

No

deg_F

Degrees Fahrenheit

Fahrenheit

No

deg_R

Rankine scale: absolute scale of thermodynamic temperature

Rankine

No

foz

U.S.

\(\mathrm{0.125\,cup}\)

fluid_oz, fluid_ounce

No

ft

International foot

\(\mathrm{12\,inch}\)

foot

No

fur

Furlong

\(\mathrm{660\,ft}\)

furlong

No

gallon

U.S.

\(\mathrm{3.7854118\,\mathcal{l}}\)

No

hp

Electrical horsepower

\(\mathrm{745.69987\,W}\)

horsepower

No

inch

International inch

\(\mathrm{2.54\,cm}\)

No

kcal

Calorie: colloquial definition of Calorie

\(\mathrm{1000\,cal}\)

Cal, Calorie, kilocal, kilocalorie

No

kip

Kilopound: force

\(\mathrm{1000\,lbf}\)

kilopound

No

kn

nautical unit of speed: 1 nmi per hour

\(\mathrm{\frac{nmi}{h}}\)

kt, knot, NMPH

No

lb

International avoirdupois pound: mass

\(\mathrm{16\,oz}\)

lbm, pound

No

lbf

Pound: force

\(\mathrm{\frac{ft\,slug}{s^{2}}}\)

No

mi

International mile

\(\mathrm{5280\,ft}\)

mile

No

mil

Thousandth of an inch

\(\mathrm{0.001\,inch}\)

thou

No

nmi

Nautical mile

\(\mathrm{1852\,m}\)

nauticalmile, NM

No

oz

International avoirdupois ounce: mass

\(\mathrm{28.349523\,g}\)

ounce

No

pint

U.S.

\(\mathrm{0.5\,quart}\)

No

psi

Pound per square inch: pressure

\(\mathrm{\frac{lbf}{inch^{2}}}\)

No

quart

U.S.

\(\mathrm{0.25\,gallon}\)

No

slug

slug: mass

\(\mathrm{32.174049\,lb}\)

No

st

International avoirdupois stone: mass

\(\mathrm{14\,lb}\)

stone

No

tbsp

U.S.

\(\mathrm{0.5\,foz}\)

tablespoon

No

ton

International avoirdupois ton: mass

\(\mathrm{2000\,lb}\)

No

tsp

U.S.

\(\mathrm{0.33333333\,tbsp}\)

teaspoon

No

yd

International yard

\(\mathrm{3\,ft}\)

yard

No

Functions

enable()

Enable Imperial units so they appear in results of find_equivalent_units and compose.

astropy.units.cds Module

This package defines units used in the CDS format, both the units defined in Centre de Données astronomiques de Strasbourg Standards for Astronomical Catalogues 2.0 format and the complete set of supported units. This format is used by VOTable up to version 1.2.

These units are not available in the top-level astropy.units namespace. To use these units, you must import the astropy.units.cds module:

>>> from astropy.units import cds
>>> q = 10. * cds.lyr  

To include them in compose and the results of find_equivalent_units, do:

>>> from astropy.units import cds
>>> cds.enable()  
Available Units

Unit

Description

Represents

Aliases

SI Prefixes

%

percent

\(\mathrm{\%}\)

No

---

dimensionless and unscaled

\(\mathrm{}\)

-

No

\h

Planck constant

\(\mathrm{6.6260701 \times 10^{-34}\,J\,s}\)

Yes

A

Ampere

\(\mathrm{A}\)

Yes

a

year

\(\mathrm{a}\)

Yes

a0

Bohr radius

\(\mathrm{5.2917721 \times 10^{-11}\,m}\)

Yes

AA

Angstrom

\(\mathrm{\mathring{A}}\)

Å, Angstrom, Angstroem

Yes

al

Light year

\(\mathrm{lyr}\)

Yes

alpha

Fine structure constant

\(\mathrm{0.0072973526\,}\)

Yes

arcmin

minute of arc

\(\mathrm{{}^{\prime}}\)

arcm

Yes

arcsec

second of arc

\(\mathrm{{}^{\prime\prime}}\)

arcs

Yes

atm

atmosphere

\(\mathrm{101325\,Pa}\)

Yes

AU

astronomical unit

\(\mathrm{AU}\)

au

Yes

bar

bar

\(\mathrm{bar}\)

Yes

barn

barn

\(\mathrm{barn}\)

Yes

bit

bit

\(\mathrm{bit}\)

Yes

byte

byte

\(\mathrm{byte}\)

Yes

C

Coulomb

\(\mathrm{C}\)

Yes

c

speed of light

\(\mathrm{2.9979246 \times 10^{8}\,\frac{m}{s}}\)

Yes

cal

calorie

\(\mathrm{4.1854\,J}\)

Yes

cd

candela

\(\mathrm{cd}\)

Yes

Crab

Crab (X-ray) flux

Yes

ct

count

\(\mathrm{ct}\)

Yes

D

Debye (dipole)

\(\mathrm{D}\)

Yes

d

Julian day

\(\mathrm{d}\)

Yes

deg

degree

\(\mathrm{{}^{\circ}}\)

°, degree

Yes

dyn

dyne

\(\mathrm{dyn}\)

Yes

e

electron charge

\(\mathrm{1.6021766 \times 10^{-19}\,C}\)

Yes

eps0

electric constant

\(\mathrm{8.8541878 \times 10^{-12}\,\frac{F}{m}}\)

Yes

erg

erg

\(\mathrm{erg}\)

Yes

eV

electron volt

\(\mathrm{eV}\)

Yes

F

Farad

\(\mathrm{F}\)

Yes

G

Gravitation constant

\(\mathrm{6.6743 \times 10^{-11}\,\frac{m^{3}}{kg\,s^{2}}}\)

Yes

g

gram

\(\mathrm{g}\)

Yes

gauss

Gauss

\(\mathrm{G}\)

Yes

geoMass

Earth mass

\(\mathrm{M_{\oplus}}\)

Mgeo

Yes

H

Henry

\(\mathrm{H}\)

Yes

h

hour

\(\mathrm{h}\)

Yes

hr

hour

\(\mathrm{h}\)

Yes

Hz

Hertz

\(\mathrm{Hz}\)

Yes

inch

inch

\(\mathrm{0.0254\,m}\)

Yes

J

Joule

\(\mathrm{J}\)

Yes

JD

Julian day

\(\mathrm{d}\)

Yes

jovMass

Jupiter mass

\(\mathrm{M_{\rm J}}\)

Mjup

Yes

Jy

Jansky

\(\mathrm{Jy}\)

Yes

K

Kelvin

\(\mathrm{K}\)

Yes

k

Boltzmann

\(\mathrm{1.380649 \times 10^{-23}\,\frac{J}{K}}\)

Yes

l

litre

\(\mathrm{\mathcal{l}}\)

Yes

lm

lumen

\(\mathrm{lm}\)

Yes

Lsun

solar luminosity

\(\mathrm{L_{\odot}}\)

solLum

Yes

lx

lux

\(\mathrm{lx}\)

Yes

lyr

Light year

\(\mathrm{lyr}\)

Yes

m

meter

\(\mathrm{m}\)

Yes

mag

magnitude

\(\mathrm{mag}\)

Yes

mas

millisecond of arc

\(\mathrm{marcsec}\)

No

me

electron mass

\(\mathrm{9.1093837 \times 10^{-31}\,kg}\)

Yes

min

minute

\(\mathrm{min}\)

Yes

MJD

Julian day

\(\mathrm{d}\)

Yes

mmHg

millimeter of mercury

\(\mathrm{133.32239\,Pa}\)

Yes

mol

mole

\(\mathrm{mol}\)

Yes

mp

proton mass

\(\mathrm{1.6726219 \times 10^{-27}\,kg}\)

Yes

Msun

solar mass

\(\mathrm{M_{\odot}}\)

solMass

Yes

mu0

magnetic constant

\(\mathrm{1.2566371 \times 10^{-6}\,\frac{N}{A^{2}}}\)

µ0

Yes

muB

Bohr magneton

\(\mathrm{9.2740101 \times 10^{-24}\,\frac{J}{T}}\)

Yes

N

Newton

\(\mathrm{N}\)

Yes

Ohm

Ohm

\(\mathrm{\Omega}\)

Yes

Pa

Pascal

\(\mathrm{Pa}\)

Yes

pc

parsec

\(\mathrm{pc}\)

Yes

ph

photon

\(\mathrm{ph}\)

Yes

pi

π

\(\mathrm{3.1415927\,}\)

Yes

pix

pixel

\(\mathrm{pix}\)

Yes

ppm

parts per million

\(\mathrm{1 \times 10^{-6}\,}\)

Yes

R

gas constant

\(\mathrm{8.3144626\,\frac{J}{K\,mol}}\)

Yes

rad

radian

\(\mathrm{rad}\)

Yes

Rgeo

Earth equatorial radius

\(\mathrm{6378100\,m}\)

Yes

Rjup

Jupiter equatorial radius

\(\mathrm{71492000\,m}\)

Yes

Rsun

solar radius

\(\mathrm{R_{\odot}}\)

solRad

Yes

Ry

Rydberg

\(\mathrm{R_{\infty}}\)

Yes

S

Siemens

\(\mathrm{S}\)

Yes

s

second

\(\mathrm{s}\)

sec

Yes

sr

steradian

\(\mathrm{sr}\)

Yes

Sun

solar unit

\(\mathrm{Sun}\)

Yes

T

Tesla

\(\mathrm{T}\)

Yes

t

metric tonne

\(\mathrm{1000\,kg}\)

Yes

u

atomic mass

\(\mathrm{1.6605391 \times 10^{-27}\,kg}\)

Yes

V

Volt

\(\mathrm{V}\)

Yes

W

Watt

\(\mathrm{W}\)

Yes

Wb

Weber

\(\mathrm{Wb}\)

Yes

yr

year

\(\mathrm{a}\)

Yes

µas

microsecond of arc

\(\mathrm{\mu arcsec}\)

No

Functions

enable()

Enable CDS units so they appear in results of find_equivalent_units and compose.

astropy.units.physical Module

Defines the physical types that correspond to different units.

Defined Physical Types

Physical type

Unit

Other physical type(s) with same unit

absement

\(\mathrm{m\,s}\)

absity

\(\mathrm{s^{2}\,m}\)

acceleration

\(\mathrm{\frac{m}{s^{2}}}\)

action

\(\mathrm{\frac{m^{2}\,kg}{s}}\)

angular momentum

amount of substance

\(\mathrm{mol}\)

angle

\(\mathrm{rad}\)

angular acceleration

\(\mathrm{\frac{rad}{s^{2}}}\)

angular frequency

\(\mathrm{\frac{rad}{s}}\)

angular speed, angular velocity

angular momentum

\(\mathrm{\frac{m^{2}\,kg}{s}}\)

action

angular speed

\(\mathrm{\frac{rad}{s}}\)

angular frequency, angular velocity

angular velocity

\(\mathrm{\frac{rad}{s}}\)

angular frequency, angular speed

area

\(\mathrm{m^{2}}\)

bandwidth

\(\mathrm{\frac{bit}{s}}\)

catalytic activity

\(\mathrm{\frac{mol}{s}}\)

chemical potential

\(\mathrm{\frac{J}{mol}}\)

column density

\(\mathrm{\frac{1}{m^{2}}}\)

compressibility

\(\mathrm{\frac{1}{Pa}}\)

crackle

\(\mathrm{\frac{m}{s^{5}}}\)

data quantity

\(\mathrm{bit}\)

diffusivity

\(\mathrm{\frac{m^{2}}{s}}\)

kinematic viscosity

dimensionless

\(\mathrm{}\)

dynamic viscosity

\(\mathrm{\frac{g}{m\,s}}\)

electrical capacitance

\(\mathrm{F}\)

electrical charge

\(\mathrm{C}\)

electrical charge (EMU)

\(\mathrm{abC}\)

electrical charge (ESU)

\(\mathrm{Fr}\)

electrical charge density

\(\mathrm{\frac{C}{m^{3}}}\)

electrical conductance

\(\mathrm{S}\)

electrical conductivity

\(\mathrm{\frac{S}{m}}\)

electrical current

\(\mathrm{A}\)

electrical current (EMU)

\(\mathrm{Bi}\)

electrical current (ESU)

\(\mathrm{statA}\)

electrical current density

\(\mathrm{\frac{A}{m^{2}}}\)

electrical dipole moment

\(\mathrm{C\,m}\)

electrical field strength

\(\mathrm{\frac{V}{m}}\)

electrical flux density

\(\mathrm{\frac{C}{m^{2}}}\)

polarization density, surface charge density

electrical impedance

\(\mathrm{\Omega}\)

electrical reactance, electrical resistance

electrical mobility

\(\mathrm{\frac{m^{2}}{V\,s}}\)

electrical potential

\(\mathrm{V}\)

electrical reactance

\(\mathrm{\Omega}\)

electrical impedance, electrical resistance

electrical resistance

\(\mathrm{\Omega}\)

electrical impedance, electrical reactance

electrical resistivity

\(\mathrm{\Omega\,m}\)

electromagnetic field strength

\(\mathrm{\frac{H}{m}}\)

permeability

electron density

\(\mathrm{\frac{e^{-}}{m^{3}}}\)

electron flux

\(\mathrm{\frac{e^{-}}{s\,m^{2}}}\)

energy

\(\mathrm{J}\)

torque, work

energy density

\(\mathrm{Pa}\)

pressure, stress

energy flux

\(\mathrm{\frac{J}{s\,m^{2}}}\)

irradiance

entropy

\(\mathrm{\frac{J}{K}}\)

heat capacity

force

\(\mathrm{N}\)

frequency

\(\mathrm{Hz}\)

frequency drift

\(\mathrm{\frac{1}{s^{2}}}\)

heat capacity

\(\mathrm{\frac{J}{K}}\)

entropy

illuminance

\(\mathrm{lx}\)

luminous emittance

impulse

\(\mathrm{\frac{kg\,m}{s}}\)

momentum

inductance

\(\mathrm{H}\)

irradiance

\(\mathrm{\frac{J}{s\,m^{2}}}\)

energy flux

jerk

\(\mathrm{\frac{m}{s^{3}}}\)

jolt

jolt

\(\mathrm{\frac{m}{s^{3}}}\)

jerk

jounce

\(\mathrm{\frac{m}{s^{4}}}\)

snap

kinematic viscosity

\(\mathrm{\frac{m^{2}}{s}}\)

diffusivity

length

\(\mathrm{m}\)

linear density

\(\mathrm{\frac{kg}{m}}\)

luminance

\(\mathrm{\frac{cd}{m^{2}}}\)

luminous efficacy

\(\mathrm{\frac{lm}{W}}\)

luminous emittance

\(\mathrm{lx}\)

illuminance

luminous flux

\(\mathrm{lm}\)

luminous intensity

\(\mathrm{cd}\)

magnetic field strength

\(\mathrm{\frac{A}{m}}\)

magnetic flux

\(\mathrm{Wb}\)

magnetic flux density

\(\mathrm{T}\)

magnetic moment

\(\mathrm{m^{2}\,A}\)

magnetic reluctance

\(\mathrm{\frac{1}{H}}\)

mass

\(\mathrm{g}\)

mass attenuation coefficient

\(\mathrm{\frac{m^{2}}{kg}}\)

opacity

mass density

\(\mathrm{\frac{kg}{m^{3}}}\)

mass flux

\(\mathrm{\frac{kg}{s\,m^{2}}}\)

momentum density

molality

\(\mathrm{\frac{mol}{kg}}\)

molar concentration

\(\mathrm{\frac{mol}{m^{3}}}\)

molar conductivity

\(\mathrm{\frac{m^{2}\,S}{mol}}\)

molar heat capacity

\(\mathrm{\frac{J}{K\,mol}}\)

molar volume

\(\mathrm{\frac{m^{3}}{mol}}\)

moment of inertia

\(\mathrm{m^{2}\,kg}\)

momentum

\(\mathrm{\frac{kg\,m}{s}}\)

impulse

momentum density

\(\mathrm{\frac{kg}{s\,m^{2}}}\)

mass flux

number density

\(\mathrm{\frac{1}{m^{3}}}\)

opacity

\(\mathrm{\frac{m^{2}}{kg}}\)

mass attenuation coefficient

particle flux

\(\mathrm{\frac{1}{s\,m^{2}}}\)

permeability

\(\mathrm{\frac{H}{m}}\)

electromagnetic field strength

permittivity

\(\mathrm{\frac{F}{m}}\)

photon flux

\(\mathrm{R}\)

photon flux density

\(\mathrm{\frac{ph}{Hz\,s\,cm^{2}}}\)

photon flux density wav

\(\mathrm{\frac{ph}{\mathring{A}\,s\,cm^{2}}}\)

plate scale

\(\mathrm{\frac{rad}{m}}\)

polarization density

\(\mathrm{\frac{C}{m^{2}}}\)

electrical flux density, surface charge density

pop

\(\mathrm{\frac{m}{s^{6}}}\)

pounce

pounce

\(\mathrm{\frac{m}{s^{6}}}\)

pop

power

\(\mathrm{W}\)

radiant flux

power density

\(\mathrm{\frac{J}{s\,m^{3}}}\)

spectral flux density wav

pressure

\(\mathrm{Pa}\)

energy density, stress

radiance

\(\mathrm{\frac{W}{sr\,m^{2}}}\)

radiant flux

\(\mathrm{W}\)

power

radiant intensity

\(\mathrm{\frac{W}{sr}}\)

reaction rate

\(\mathrm{\frac{mol}{s\,m^{3}}}\)

snap

\(\mathrm{\frac{m}{s^{4}}}\)

jounce

solid angle

\(\mathrm{sr}\)

specific energy

\(\mathrm{\frac{J}{kg}}\)

specific entropy

\(\mathrm{\frac{J}{K\,kg}}\)

specific heat capacity

specific heat capacity

\(\mathrm{\frac{J}{K\,kg}}\)

specific entropy

specific volume

\(\mathrm{\frac{m^{3}}{kg}}\)

spectral flux density

\(\mathrm{Jy}\)

spectral flux density wav

\(\mathrm{\frac{J}{s\,m^{3}}}\)

power density

speed

\(\mathrm{\frac{m}{s}}\)

velocity

stress

\(\mathrm{Pa}\)

energy density, pressure

surface charge density

\(\mathrm{\frac{C}{m^{2}}}\)

electrical flux density, polarization density

surface mass density

\(\mathrm{\frac{kg}{m^{2}}}\)

surface tension

\(\mathrm{\frac{m^{2}\,W}{Hz}}\)

temperature

\(\mathrm{K}\)

temperature gradient

\(\mathrm{\frac{K}{m}}\)

thermal conductance

\(\mathrm{\frac{W}{K}}\)

thermal conductivity

\(\mathrm{\frac{W}{K\,m}}\)

thermal resistance

\(\mathrm{\frac{K}{W}}\)

thermal resistivity

\(\mathrm{\frac{K\,m}{W}}\)

time

\(\mathrm{s}\)

torque

\(\mathrm{J}\)

energy, work

velocity

\(\mathrm{\frac{m}{s}}\)

speed

volume

\(\mathrm{m^{3}}\)

volumetric flow rate

\(\mathrm{\frac{m^{3}}{s}}\)

volumetric rate

\(\mathrm{\frac{1}{s\,m^{3}}}\)

wavenumber

\(\mathrm{\frac{1}{m}}\)

work

\(\mathrm{J}\)

energy, torque

yank

\(\mathrm{\frac{N}{s}}\)

Functions

def_physical_type(unit, name)

Add a mapping between a unit and the corresponding physical type(s).

get_physical_type(obj)

Return the physical type that corresponds to a unit (or another physical type representation).

Classes

PhysicalType(unit, physical_types)

Represents the physical type(s) that are dimensionally compatible with a set of units.

Class Inheritance Diagram

Inheritance diagram of astropy.units.physical.PhysicalType

astropy.units.equivalencies Module

A set of standard astronomical equivalencies.

Functions

parallax()

Returns a list of equivalence pairs that handle the conversion between parallax angle and distance.

spectral()

Returns a list of equivalence pairs that handle spectral wavelength, wave number, frequency, and energy equivalencies.

spectral_density(wav[, factor])

Returns a list of equivalence pairs that handle spectral density with regard to wavelength and frequency.

doppler_radio(rest)

Return the equivalency pairs for the radio convention for velocity.

doppler_optical(rest)

Return the equivalency pairs for the optical convention for velocity.

doppler_relativistic(rest)

Return the equivalency pairs for the relativistic convention for velocity.

doppler_redshift()

Returns the equivalence between Doppler redshift (unitless) and radial velocity.

mass_energy()

Returns a list of equivalence pairs that handle the conversion between mass and energy.

brightness_temperature(frequency[, beam_area])

Defines the conversion between Jy/sr and "brightness temperature", \(T_B\), in Kelvins.

thermodynamic_temperature(frequency[, T_cmb])

Defines the conversion between Jy/sr and "thermodynamic temperature", \(T_{CMB}\), in Kelvins.

beam_angular_area(beam_area)

Convert between the beam unit, which is commonly used to express the area of a radio telescope resolution element, and an area on the sky.

dimensionless_angles()

Allow angles to be equivalent to dimensionless (with 1 rad = 1 m/m = 1).

logarithmic()

Allow logarithmic units to be converted to dimensionless fractions

temperature()

Convert between Kelvin, Celsius, Rankine and Fahrenheit here because Unit and CompositeUnit cannot do addition or subtraction properly.

temperature_energy()

Convert between Kelvin and keV(eV) to an equivalent amount.

molar_mass_amu()

Returns the equivalence between amu and molar mass.

pixel_scale(pixscale)

Convert between pixel distances (in units of pix) and other units, given a particular pixscale.

plate_scale(platescale)

Convert between lengths (to be interpreted as lengths in the focal plane) and angular units with a specified platescale.

Classes

Equivalency(equiv_list[, name, kwargs])

A container for a units equivalency.

Class Inheritance Diagram

Inheritance diagram of astropy.units.equivalencies.Equivalency

astropy.units.function.core Module

Function Units and Quantities.

Classes

FunctionUnitBase([physical_unit, function_unit])

Abstract base class for function units.

FunctionQuantity(value[, unit, dtype, copy, ...])

A representation of a (scaled) function of a number with a unit.

Class Inheritance Diagram

Inheritance diagram of astropy.units.function.core.FunctionUnitBase, astropy.units.function.core.FunctionQuantity

astropy.units.function.logarithmic Module

Classes

LogUnit([physical_unit, function_unit])

Logarithmic unit containing a physical one

MagUnit([physical_unit, function_unit])

Logarithmic physical units expressed in magnitudes

DexUnit([physical_unit, function_unit])

Logarithmic physical units expressed in magnitudes

DecibelUnit([physical_unit, function_unit])

Logarithmic physical units expressed in dB

LogQuantity(value[, unit, dtype, copy, ...])

A representation of a (scaled) logarithm of a number with a unit

Magnitude(value[, unit, dtype, copy, order, ...])

Decibel(value[, unit, dtype, copy, order, ...])

Dex(value[, unit, dtype, copy, order, ...])

Variables

STmag

ST magnitude: STmag=-21.1 corresponds to 1 erg/s/cm2/A

ABmag

AB magnitude: ABmag=-48.6 corresponds to 1 erg/s/cm2/Hz

M_bol

Absolute bolometric magnitude: M_bol=0 corresponds to L_bol0=3.0128e+28 J / s

m_bol

Apparent bolometric magnitude: m_bol=0 corresponds to f_bol0=2.51802e-08 kg / s3

Class Inheritance Diagram

Inheritance diagram of astropy.units.function.logarithmic.LogUnit, astropy.units.function.logarithmic.MagUnit, astropy.units.function.logarithmic.DexUnit, astropy.units.function.logarithmic.DecibelUnit, astropy.units.function.logarithmic.LogQuantity, astropy.units.function.logarithmic.Magnitude, astropy.units.function.logarithmic.Decibel, astropy.units.function.logarithmic.Dex

astropy.units.deprecated Module

This package defines deprecated units.

These units are not available in the top-level astropy.units namespace. To use these units, you must import the astropy.units.deprecated module:

>>> from astropy.units import deprecated
>>> q = 10. * deprecated.emu  

To include them in compose and the results of find_equivalent_units, do:

>>> from astropy.units import deprecated
>>> deprecated.enable()  
Available Units

Unit

Description

Represents

Aliases

SI Prefixes

emu

Biot: CGS (EMU) unit of current

\(\mathrm{Bi}\)

No

Prefixes for earthMass

Earth mass prefixes

\(\mathrm{5.9721679 \times 10^{24}\,kg}\)

M_earth, Mearth

Only

Prefixes for earthRad

Earth radius prefixes

\(\mathrm{6378100\,m}\)

R_earth, Rearth

Only

Prefixes for jupiterMass

Jupiter mass prefixes

\(\mathrm{1.8981246 \times 10^{27}\,kg}\)

M_jup, Mjup, M_jupiter, Mjupiter

Only

Prefixes for jupiterRad

Jupiter radius prefixes

\(\mathrm{71492000\,m}\)

R_jup, Rjup, R_jupiter, Rjupiter

Only

Functions

enable()

Enable deprecated units so they appear in results of find_equivalent_units and compose.

astropy.units.required_by_vounit Module

This package defines SI prefixed units that are required by the VOUnit standard but that are rarely used in practice and liable to lead to confusion (such as msolMass for milli-solar mass). They are in a separate module from astropy.units.deprecated because they need to be enabled by default for astropy.units to parse compliant VOUnit strings. As a result, e.g., Unit('msolMass') will just work, but to access the unit directly, use astropy.units.required_by_vounit.msolMass instead of the more typical idiom possible for the non-prefixed unit, astropy.units.solMass.

Available Units

Unit

Description

Represents

Aliases

SI Prefixes

Prefixes for solLum

Solar luminance prefixes

\(\mathrm{3.828 \times 10^{26}\,W}\)

L_sun, Lsun

Only

Prefixes for solMass

Solar mass prefixes

\(\mathrm{1.9884099 \times 10^{30}\,kg}\)

M_sun, Msun

Only

Prefixes for solRad

Solar radius prefixes

\(\mathrm{6.957 \times 10^{8}\,m}\)

R_sun, Rsun

Only