astropy:docs

UnitBase

class astropy.units.UnitBase[source] [edit on github]

Bases: object

Abstract base class for units.

Most of the arithmetic operations on units are defined in this base class.

Should not be instantiated by users directly.

Attributes Summary

aliases Returns the alias (long) names for this unit.
bases Return the bases of the unit.
cgs Returns a copy of the current Unit instance with CGS units.
name Returns the canonical (short) name associated with this unit.
names Returns all of the names associated with this unit.
physical_type Return the physical type on the unit.
powers Return the powers of the unit.
scale Return the scale of the unit.
si Returns a copy of the current Unit instance in SI units.

Methods Summary

compose([equivalencies, units, max_depth, ...]) Return the simplest possible composite unit(s) that represent the given unit.
decompose([bases]) Return a unit object composed of only irreducible units.
find_equivalent_units([equivalencies, ...]) Return a list of all the units that are the same type as self.
get_converter(*args, **kwargs)

Deprecated since version 1.0.

in_units(other[, value, equivalencies]) Alias for to for backward compatibility with pynbody.
is_equivalent(other[, equivalencies]) Returns True if this unit is equivalent to other.
is_unity() Returns True if the unit is unscaled and dimensionless.
to(other[, value, equivalencies]) Return the converted values in the specified unit.
to_string([format]) Output the unit in the given format as a string.
to_system(system) Converts this unit into ones belonging to the given system.

Attributes Documentation

aliases

Returns the alias (long) names for this unit.

bases

Return the bases of the unit.

cgs

Returns a copy of the current Unit instance with CGS units.

name

Returns the canonical (short) name associated with this unit.

names

Returns all of the names associated with this unit.

physical_type

Return the physical type on the unit.

Examples

>>> from astropy import units as u
>>> print(u.m.physical_type)
length
powers

Return the powers of the unit.

scale

Return the scale of the unit.

si

Returns a copy of the current Unit instance in SI units.

Methods Documentation

compose(equivalencies=[], units=None, max_depth=2, include_prefix_units=False)[source] [edit on github]

Return the simplest possible composite unit(s) that represent the given unit. Since there may be multiple equally simple compositions of the unit, a list of units is always returned.

Parameters:

equivalencies : list of equivalence pairs, optional

A list of equivalence pairs to also list. See Equivalencies. This list is in addition to possible global defaults set by, e.g., set_enabled_equivalencies. Use None to turn off all equivalencies.

units : set of units to compose to, optional

If not provided, any known units may be used to compose into. Otherwise, units is a dict, module or sequence containing the units to compose into.

max_depth : int, optional

The maximum recursion depth to use when composing into composite units.

include_prefix_units : bool, optional

When True, include prefixed units in the result. Default is False.

Returns:

units : list of CompositeUnit

A list of candidate compositions. These will all be equally simple, but it may not be possible to automatically determine which of the candidates are better.

decompose(bases=set([]))[source] [edit on github]

Return a unit object composed of only irreducible units.

Parameters:

bases : sequence of UnitBase, optional

The bases to decompose into. When not provided, decomposes down to any irreducible units. When provided, the decomposed result will only contain the given units. This will raises a UnitsError if it’s not possible to do so.

Returns:

unit : CompositeUnit object

New object containing only irreducible unit objects.

find_equivalent_units(equivalencies=[], units=None, include_prefix_units=False)[source] [edit on github]

Return a list of all the units that are the same type as self.

Parameters:

equivalencies : list of equivalence pairs, optional

A list of equivalence pairs to also list. See Equivalencies. Any list given, including an empty one, supercedes global defaults that may be in effect (as set by set_enabled_equivalencies)

units : set of units to search in, optional

If not provided, all defined units will be searched for equivalencies. Otherwise, may be a dict, module or sequence containing the units to search for equivalencies.

include_prefix_units : bool, optional

When True, include prefixed units in the result. Default is False.

Returns:

units : list of UnitBase

A list of unit objects that match u. A subclass of list (EquivalentUnitsList) is returned that pretty-prints the list of units when output.

get_converter(*args, **kwargs)[source] [edit on github]

Deprecated since version 1.0: The get_converter function is deprecated and may be removed in a future version.

Return the conversion function to convert values from self to the specified unit.

Parameters:

other : unit object or string

The unit to convert to.

equivalencies : list of equivalence pairs, optional

A list of equivalence pairs to try if the units are not directly convertible. See Equivalencies. This list is in addition to possible global defaults set by, e.g., set_enabled_equivalencies. Use None to turn off all equivalencies.

Returns:

func : callable

A callable that normally expects a single argument that is a scalar value or an array of values (or anything that may be converted to an array).

Raises:

UnitsError

If units are inconsistent

in_units(other, value=1.0, equivalencies=[])[source] [edit on github]

Alias for to for backward compatibility with pynbody.

is_equivalent(other, equivalencies=[])[source] [edit on github]

Returns True if this unit is equivalent to other.

Parameters:

other : unit object or string or tuple

The unit to convert to. If a tuple of units is specified, this method returns true if the unit matches any of those in the tuple.

equivalencies : list of equivalence pairs, optional

A list of equivalence pairs to try if the units are not directly convertible. See Equivalencies. This list is in addition to possible global defaults set by, e.g., set_enabled_equivalencies. Use None to turn off all equivalencies.

Returns:

bool

is_unity()[source] [edit on github]

Returns True if the unit is unscaled and dimensionless.

to(other, value=1.0, equivalencies=[])[source] [edit on github]

Return the converted values in the specified unit.

Parameters:

other : unit object or string

The unit to convert to.

value : scalar int or float, or sequence convertable to array, optional

Value(s) in the current unit to be converted to the specified unit. If not provided, defaults to 1.0

equivalencies : list of equivalence pairs, optional

A list of equivalence pairs to try if the units are not directly convertible. See Equivalencies. This list is in addition to possible global defaults set by, e.g., set_enabled_equivalencies. Use None to turn off all equivalencies.

Returns:

values : scalar or array

Converted value(s). Input value sequences are returned as numpy arrays.

Raises:

UnitsError

If units are inconsistent

to_string(format=<class 'astropy.units.format.generic.Generic'>)[source] [edit on github]

Output the unit in the given format as a string.

Parameters:

format : astropy.units.format.Base instance or str

The name of a format or a formatter object. If not provided, defaults to the generic format.

to_system(system)[source] [edit on github]

Converts this unit into ones belonging to the given system. Since more than one result may be possible, a list is always returned.

Parameters:

system : module

The module that defines the unit system. Commonly used ones include astropy.units.si and astropy.units.cgs.

To use your own module it must contain unit objects and a sequence member named bases containing the base units of the system.

Returns:

units : list of CompositeUnit

The list is ranked so that units containing only the base units of that system will appear first.

Page Contents