Package ivs :: Package units :: Module conversions
[hide private]
[frames] | no frames]

Module conversions

source code

Convert one unit (and uncertainty) to another.

Contents:

  1. The Python module: basic usage of the Python module
  2. The Terminal tool: basic usage of the terminal tool
  3. Fundamental constants and base unit systems
  4. Calculating with units
  5. Dealing with and interpreting unit strings

Much of the unit conversions has been tested using Appendix B of http://physics.nist.gov/cuu/pdf/sp811.pdf, and examples within, which is the international document describing the SI standard. As a matter of fact, this module should be fully compatible with the SI unit conventions, except for the notation of units (brackets are not allowed). This module extends the SI unit conventions to be more flexible and intuitive, and also allows to ditch the SI unit convention alltogether.

Some of the many possibilities include (see convert for an extensive set of examples):

  1. Conversions between equal-type units: meter to nano-lightyears, erg/s to W, cy/d to muHz, but also erg/s/cm2/A to W/m2/mum, sr to deg2, etc...
  2. Conversions between unequal-type units: angstrom to km/s via the speed of light, F(lambda) to F(nu), F(nu) to lambdaF(lambda)/sr, meter to cycles/arcsec (interferometry), etc...
  3. Nonlinear conversions: vegamag to erg/s/cm2/A or Jy, Celcius to Fahrenheit or Kelvin, calender date to (any kind of modified) Julian Day, ...
  4. Conversions of magnitude to flux amplitudes via 'Amag' and 'ppt' or 'ampl'
  5. Conversions of magnitude colors to flux ratios via 'mag_color' and 'flux_ratio'
  6. Coordinate transformations (equatorial to galactic etc.)
  7. Logarithmic conversions, e.g. from logTeff to Teff via '[K]' and 'K'
  8. Inclusion of uncertainties, both in input values and/or reference values when converting between unequal-type units, automatically recognised when giving two positional argument (value, error) instead of one (value).
  9. Currency exchange rates. First, you need to call set_exchange_rates for this to work, since the latest currency definitions and rates need to queried from the European Central Bank (automatically done when using the terminal tool).
  10. Computations with units.

Warning 1: frequency units are technically given in cycles per time (cy). This means that if you want to convert e.g. muHz to d-1 (or 1/d), you need to ask for cy/d. There is a general ambiguity concerning the unit 'cycles': it is not an official SI unit, but is needed to make the distinction between e.g. rad/s and Hz. In true SI-style, the official "unit" of rad/s is actually just s-1, since radians are not really a unit. This gives room for confusion! To make everything even more confusing, there is another unit which is equal to the reciprocal second, namely the Becquerel. This is used for stochastic or non-recurrent phenomena. Basically, the problem is:

   rad/s == 1/s
   1/s   == Hz

but:

   rad/s != Hz

If you have any doubts, set the logger to display debug messages (e.g. logger.setLevel(10), see ivs.aux.loggers). It will sometimes tell you if certain assumptions are made.

Warning 2: there is some ambiguity between units. For example, as can be interpreted as 'arcsecond', but also as 'attosecond'. In case of ambiguity, the basic units will be preferred (in this case 'arcsecond'), instead of the one with prefixes. This is a list of non-exhaustive ambiguities (the preffered one in italic):

Sometimes the case-sensitivity of the conversions comes in handy. There is no ambiguity between:

Warning 3: the unit name of angstrom is AA, ampere is A.

Warning 4: Most of the imperial units are UK/Canada. If you need US, prefix the unit with US: E.g. The gallon (gal) is the international imperial gallon, USgal is the US gallon.

Note 1: Photometric passbands are given as a string "SYSTEM.FILTER". For a list of defined systems and filters, see ivs.sed.filters.

Section 1. The Python module

The main function convert (see link for a full list of examples) does all the work and is called via

>>> result = convert('km','m',1.)

or when the error is known

>>> result,e_result = convert('km','m',1.,0.1)

Be careful when mixing nonlinear conversions (e.g. magnitude to flux) with linear conversions (e.g. Jy to W/m2/m).

When your favorite conversion is not implemented, there are five places where you can add information:

  1. _scalings: if your favorite prefix (e.g. Tera, nano...) is not available
  2. _aliases: if your unit is available but not under the name you are used to.
  3. _factors: if your unit is not available.
  4. _switch: if your units are available, but conversions from one to another is not straightforward and extra infromation is needed (e.g. to go from angstrom to km/s in a spectrum, a reference wavelength is needed).
  5. _convention: if your favorite base unit system is not defined (SI,cgs...).

If you need to add a linear factor, just give the factor in SI units, and the SI base units it consists of (see _factors for examples). If you need to add a nonlinear factor, you need to give a class definition (see the examples).

Section 2. The Terminal tool

For help and a list of all defined units and abbreviations, do:

   $:> python conversions.py --help
   =====================================           | =====================================           
   =   Units of absorbed dose          =           | =   Units of acceleration           =           
   =====================================           | =====================================           
   Gy              = gray                          | Gal             = Gal                           
   =====================================           | =====================================           
   =   Units of angle                  =           | =   Units of area                   =           
   =====================================           | =====================================           
   am              = arcminute                     | a               = are                           
   as              = arcsecond                     | ac              = acre (international)          
   cy              = cycle                         | b               = barn                          
   deg             = degree                        | =====================================           
   rad             = radian                        | =   Units of coordinate             =           
   rpm             = revolutions per minute        | =====================================           
   sr              = sterradian                    | complex_coord   = <own unit>                    
   =====================================           | deg_coord       = degrees                       
   =   Units of catalytic activity     =           | ecliptic        = ecliptic                      
   =====================================           | equatorial      = equatorial                    
   kat             = katal                         | galactic        = galactic                      
   =====================================           | rad_coord       = radians                       
   =   Units of currency               =           | =====================================           
   =====================================           | =   Units of dose equivalent        =           
   EUR             = EURO                          | =====================================           
   =====================================           | Sv              = sievert                       
   =   Units of dynamic viscosity      =           | rem             = rem                           
   =====================================           | =====================================           
   P               = poise                         | =   Units of electric capacitance   =           
   =====================================           | =====================================           
   =   Units of electric charge        =           | F               = Farad                         
   =====================================           | =====================================           
   C               = Coulomb                       | =   Units of electric conductance   =           
   =====================================           | =====================================           
   =   Units of electric current       =           | S               = Siemens                       
   =====================================           | =====================================           
   A               = Ampere                        | =   Units of electric potential difference   =  
   Bi              = biot                          | =====================================           
   Gi              = Ampere                        | V               = Volt                          
   =====================================           | =====================================           
   =   Units of electric resistance    =           | =   Units of energy                 =           
   =====================================           | =====================================           
   O               = Ohm                           | Cal             = large calorie (international table)
   =====================================           | J               = Joule                         
   =   Units of energy/power           =           | cal             = small calorie (international table)
   =====================================           | eV              = electron volt                 
   Lsol            = Solar luminosity              | erg             = ergon                         
   =====================================           | foe             = (ten to the) fifty one ergs   
   =   Units of flux density           =           | =====================================           
   =====================================           | =   Units of flux                   =           
   Jy              = Jansky                        | =====================================           
   =====================================           | ABmag           = AB magnitude                  
   =   Units of frequency              =           | Amag            = amplitude in magnitude        
   =====================================           | STmag           = ST magnitude                  
   hz              = Hertz                         | ampl            = fractional amplitude          
   =====================================           | flux_ratio      = flux ratio                    
   =   Units of inductance             =           | mag             = magnitude                     
   =====================================           | mag_color       = color                         
   H               = Henry                         | pph             = amplitude in parts per hundred
   =====================================           | ppm             = amplitude in parts per million
   =   Units of length                 =           | ppt             = amplitude in parts per thousand
   =====================================           | vegamag         = Vega magnitude                
   AA              = angstrom                      | =====================================           
   AU              = astronomical unit             | =   Units of force                  =           
   Rearth          = Earth radius                  | =====================================           
   Rjup            = Jupiter radius                | N               = Newton                        
   Rsol            = Solar radius                  | dyn             = dyne                          
   USft            = foot (US)                     | =====================================           
   USmi            = mile (US)                     | =   Units of illuminance            =           
   a0              = Bohr radius                   | =====================================           
   bs              = beard second                  | lx              = lux                           
   ch              = chain                         | ph              = phot                          
   ell             = ell                           | sb              = stilb                         
   fathom          = fathom                        | =====================================           
   ft              = foot (international)          | =   Units of kynamatic viscosity    =           
   fur             = furlong                       | =====================================           
   in              = inch (international)          | St              = stokes                        
   ly              = light year                    | =====================================           
   m               = meter                         | =   Units of luminous flux          =           
   mi              = mile (international)          | =====================================           
   nami            = nautical mile                 | lm              = lumen                         
   pc              = parsec                        | =====================================           
   perch           = pole                          | =   Units of magnetic field strength   =        
   pole            = perch                         | =====================================           
   potrzebie       = potrzebie                     | G               = Gauss                         
   rd              = rod                           | T               = Tesla                         
   smoot           = smooth                        | =====================================           
   yd              = yard (international)          | =   Units of magnetizing field      =           
   =====================================           | =====================================           
   =   Units of m/s                    =           | Oe              = Oersted                       
   =====================================           | =====================================           
   cc              = Speed of light                | =   Units of power                  =           
   =====================================           | =====================================           
   =   Units of magnetic flux          =           | W               = Watt                          
   =====================================           | dp              = Donkeypower                   
   Mx              = Maxwell                       | hp              = Horsepower                    
   Wb              = Weber                         | =====================================           
   =====================================           | =   Units of roentgen               =           
   =   Units of mass                   =           | =====================================           
   =====================================           | R               = Roentgen                      
   Mearth          = Earth mass                    | =====================================           
   Mjup            = Jupiter mass                  | =   Units of temperature            =           
   Mlun            = Lunar mass                    | =====================================           
   Msol            = Solar mass                    | Cel             = Celcius                       
   carat           = carat                         | Far             = Fahrenheit                    
   firkin          = firkin                        | K               = Kelvin                        
   g               = gram                          | Tsol            = Solar temperature             
   gr              = gram                          | =====================================           
   lb              = pound                         | =   Units of velocity               =           
   mol             = molar mass                    | =====================================           
   ounce           = ounce                         | knot            = nautical mile per hour        
   st              = stone                         | mph             = miles per hour                
   ton             = gram                          | 
   u               = atomic mass                   | 
   =====================================           | 
   =   Units of pressure               =           | 
   =====================================           | 
   Pa              = Pascal                        | 
   at              = atmosphere (technical)        | 
   atm             = atmosphere (standard)         | 
   ba              = barye                         | 
   bar             = baros                         | 
   mmHg            = millimeter of mercury         | 
   psi             = pound per square inch         | 
   torr            = Torricelli                    | 
   =====================================           | 
   =   Units of second                 =           | 
   =====================================           | 
   fortnight       = fortnight                     | 
   =====================================           | 
   =   Units of time                   =           | 
   =====================================           | 
   Bq              = Becquerel                     | 
   CD              = calender day                  | 
   Ci              = Curie                         | 
   JD              = Julian day                    | 
   MJD             = modified Julian day           | 
   cr              = century                       | 
   d               = day                           | 
   h               = hour                          | 
   j               = jiffy                         | 
   min             = minute                        | 
   mo              = month                         | 
   s               = second                        | 
   sidereal        = sidereal day                  | 
   wk              = week                          | 
   yr              = year                          | 
   =====================================           | 
   =   Units of volume                 =           | 
   =====================================           | 
   USgal           = gallon (US)                   | 
   USgi            = gill (Canadian and UK imperial)| 
   bbl             = barrel                        | 
   bu              = bushel                        | 
   gal             = gallon (Canadian and UK imperial)| 
   gi              = gill (Canadian and UK imperial)| 
   l               = liter                         | 
   ngogn           = 1000 cubic potrzebies         | 


   Usage: conversions.py --from=<unit> --to=<unit> [options] value [error]

   Options:
   -h, --help            show this help message and exit
   --from=_FROM          units to convert from
   --to=_TO              units to convert to

   Extra quantities when changing base units (e.g. Flambda to Fnu):
       -w WAVE, --wave=WAVE
                           wavelength with units (e.g. used to convert Flambda to
                           Fnu)
       -f FREQ, --freq=FREQ
                           frequency (e.g. used to convert Flambda to Fnu)
       -p PHOTBAND, --photband=PHOTBAND
                           photometric passband

To convert from one unit to another, do:

   $:> python conversions.py --from=nRsol/h --to cm/s=1.2345
   1.2345 nRsol/h    =    0.0238501 cm/s

In fact, the to parameter is optional (despite it not being a positional argument, from is not optional). The script will assume you want to convert to SI:

   $:> python conversions.py --from=nRsol/h 1.2345
   1.2345 nRsol/h    =    0.000238501 m1 s-1

It is also possible to compute with uncertainties, and you can give extra keyword arguments if extra information is needed for conversions:

   $:> python conversions.py --from=mag --to=erg/s/cm2/AA --photband=GENEVA.U 7.84 0.02
   7.84 +/- 0.02 mag    =    4.12191e-12 +/- 7.59283e-14 erg/s/cm2/AA

If you want to do coordinate transformations, e.g. from fractional radians to degrees/arcminutes/arcseconds, you can do:

   $:> python conversions.py --from=rad_coord --to=equatorial 5.412303,0.123
     (5.412303, 0.123) rad_coord    =    20:40:24.51,7:02:50.6 equ

Section 3. Fundamental constants and base unit systems

Although fundamental constants are supposed to be constants, there are two reasons why one might to change their value:

  1. You want to work in a different base unit system (e.g. cgs instead of SI)
  2. You disagree with some of the literature values and want to use your own.

Section 3.1. Changing base unit system

A solution to the first option might be to define all constants in SI, and define them again in cgs, adding a postfix _cgs to the constants' name. This is done in the module ivs.units.constants to give the user access to values of constants in common-used systems. However, this way of working does not offer any flexibility, and one has to redefine all the values manually for any given convention. Also, you may want to work in cgs by default, which means that typing the postfix every time is cumbersome. For this purpose, you can redefine all constants in a different base system with a simple command:

>>> print constants.Msol,constants.GG
1.988547e+30 6.67384e-11
>>> set_convention(units='cgs')
('SI', 'standard', 'rad')
>>> print constants.Msol,constants.GG
1.988547e+33 6.67384e-08

Remark that the units are also changed accordingly:

>>> print constants.Msol_units
g1

You can go crazy with this:

>>> set_convention(units='imperial')
('cgs', 'standard', 'rad')
>>> print(constants.GG,constants.GG_units)
(3.9594319112470405e-11, 'yd3 lb-1 s-2')

Resetting to the default SI can be done by calling set_convention without any arguments, or simply

>>> set_convention(units='SI')
('imperial', 'standard', 'rad')

The function set_convention returns the current settings, so you can remember the old settings and make a temporary switch:

>>> old_settings = set_convention(units='cgs')
>>> set_convention(*old_settings)
('cgs', 'standard', 'rad')

Warning: Changing the value of the constants affects all modules, where an import statement from ivs.units import constants is made. It will not change the values in modules where the import is done via from ivs.units.constants import *. If you want to get the value of a fundamental constant regardless of the preference of base system, you might want to call

>>> Msol = get_constant('Msol','cgs')

To change the behaviour of the interpretation of radians and cycle, you can specify the keyword frequency. See set_convention for more information.

Section 3.2. Changing the values of the fundamental constants

If you disagree with some of the literature values, you can also redefine the values of the fundamental constants. For example, to use the values defined in the stellar evolution code MESA, simply do

>>> set_convention(units='cgs',values='mesa')
('SI', 'standard', 'rad')

You can query for info on the current convention without changing it:

>>> get_convention()
('cgs', 'mesa', 'rad')

But for the sake of the examples, we'll switch back to the default SI...

>>> reset()

Section 4. Calculating with units

There exists a class Unit, which allows you to calculate with values with units (optionally including uncertainties). The class Unit offers high flexibility when it comes to initialization, so that a user can create a Unit according to his/her own preferences. See Unit.__init__ for all the options.

The purpose of existence of the class Unit is calculation: you can simply create a couple of units, and multiply, add, divide etc.. with Python's standard operators ('*','+','/',...). The standard rules apply: you can multiply and divide basically any combination of parameters, but adding and subtracting requires Units to have the same dimensions!

Aside from the basic operators, a Unit also understands the most commonly used operations from numpy (np.sin, np.cos, np.sqrt...). If a function is not implemented, you can easily add it yourself following the existing examples in the code. A Unit does not understand the functions from the math module, only the numpy implementation!

Better yet, you can also put Units in a numpy array. Just make a list of units and make it an array as you would do with normal floats (no need to specify the dtype). You can safely mix units with and without uncertainties, and mix units with different dimensions. You have access to most of numpy functions such as mean, sum etc. Of course, when you mix different dimensions, calling sum will result in an error!

Example 1: Calculate the inclination of star, given a measured vsini=150+/-11 km/s, a rotation period of 2+/-1 days and a radius of 12 solar radii (no errors). First we give in our data:

>>> vsini = Unit((150.,11),'km/s')
>>> P = Unit((2.,1.),'d')
>>> R = Unit(12.,'Rsol')

Now we can calculate the sini:

>>> sini =  vsini * P / (2*np.pi*R)

And take the arcsine to recover the inclination angle in radians. Because we are working with the Unit class, we can also immediately retrieve it in degrees:

>>> print np.arcsin(sini)
0.517004704077+/-0.287337185083 rad
>>> print np.arcsin(sini).convert('deg')
29.622187532+/-16.4632080024 deg

Example 2: The following is an exam question on the Biophysics exam (1st year Bachelor) about error propagation.

Question: Suppose there is a party to celebrate the end of the Biophysics exam. You want to invite 4 persons, and you want to know if 1 liter of champagne is enough to fill 5 glasses. The glasses are cylinders with a circumference of 15.5+/-0.5cm, and a height of 10.0+/-0.5cm. Calculate the volume of one glass and its uncertainty. Can you pour 5 glasses of champagne from the 1 liter bottle?

Answer:

>>> r = Unit( (15.5/(2*np.pi), 0.5/(2*np.pi)), 'cm')
>>> h = Unit( (10.0,0.5), 'cm')
>>> V = np.pi*r**2*h
>>> print(V)
0.000191184875389+/-1.56051027314e-05 m3
>>> print((5*V).convert('liter'))
0.955924376946+/-0.0780255136569 liter

It is not sufficient within about 1 sigma.

Section 5, Dealing with and interpreting string units

Suppose you want to make automated plots with a human readable name and a LaTeX string of the unit. Then you can use get_type and unit2texlabel:

>>> print(unit2texlabel('erg/s/mum2/AA',full=True))
Flux Density [erg  s$^{-1}$ $\mu$m$^{-2}$ $\AA$$^{-1}$]
>>> print(unit2texlabel('P',full=True))
Dynamic Viscosity [P]

or

>>> print(unit2texlabel('erg/s/mum2/AA'))
erg  s$^{-1}$ $\mu$m$^{-2}$ $\AA$$^{-1}$
Classes [hide private]
    Nonlinear change-of-base functions
  NonLinearConverter
Base class for nonlinear conversions
  Fahrenheit
Convert Fahrenheit to Kelvin and back
  Celcius
Convert Celcius to Kelvin and back
  AmplMag
Convert a Vega magnitude to W/m2/m (Flambda) and back
  VegaMag
Convert a Vega magnitude to W/m2/m (Flambda) and back
  ABMag
Convert an AB magnitude to W/m2/Hz (Fnu) and back
  STMag
Convert an ST magnitude to W/m2/m (Flambda) and back
  Color
Convert a color to a flux ratio and back
  DecibelSPL
Convert a Decibel to intensity W/m2 and back.
  JulianDay
Convert a calender date to Julian date and back
  ModJulianDay
Convert a Modified Julian Day to Julian Day and back
  GalCoords
Convert Galactic coords to complex coords and back
  EquCoords
Convert Equatorial coords to complex coords and back
  EclCoords
Convert Ecliptic coords to complex coords and back
  DegCoords
Convert Complex coords to degrees coordinates and back
  RadCoords
Convert Complex coords to degrees coordinates and back
    Computations with units
  Unit
Class to calculate with numbers (and uncertainties) containing units.
Functions [hide private]
    Main functions
float
convert(_from, _to, *args, **kwargs)
Convert one unit to another.
source code
array
nconvert(_froms, _tos, *args, **kwargs)
Convert a list/array/tuple of values with different units to other units.
source code
str
change_convention(to_, units, origin=None)
Change units from one convention to another.
source code
str
set_convention(units='SI', values='standard', frequency='rad')
Consistently change the values of the fundamental constants to the paradigm of other system units or programs.
source code
 
reset()
Resets all values and conventions to the original values.
source code
str,str
get_convention()
Returns convention and name of set of values.
source code
float
get_constant(constant_name, units='SI', value='standard')
Convenience function to retrieve the value of a constant in a particular system.
source code
dict
get_constants(units='SI', values='standard')
Convenience function to retrieve the value of all constants in a particular system.
source code
    Conversions basics and helper functions
str
solve_aliases(unit)
Resolve simple aliases in a unit's name.
source code
(float,str,int)
components(unit)
Decompose a unit into: factor, SI base units, power.
source code
(float,str)
breakdown(unit)
Decompose a unit into SI base units containing powers.
source code
str
compress(unit, ignore_factor=False)
Compress basic units to more human-readable type.
source code
list,list
split_units_powers(unit)
Split unit in a list of units and a list of powers.
source code
bool
is_basic_unit(unit, type)
Check if a unit is represents one of the basic quantities (length, mass...)
source code
bool
is_type(unit, type)
Check if a unit is of a certain type (mass, length, force...)
source code
str
get_type(unit)
Return human readable name of the unit
source code
float
round_arbitrary(x, base=5)
Round to an arbitrary base.
source code
str
unit2texlabel(unit, full=False)
Convert a unit string to a nicely formatted TeX label.
source code
str
get_help()
Return a string with a list and explanation of all defined units.
source code
 
units2sphinx() source code
 
derive_wrapper(fctn) source code
    Linear change-of-base conversions
float
distance2velocity(arg, **kwargs)
Switch from distance to velocity via a reference wavelength.
source code
float
velocity2distance(arg, **kwargs)
Switch from velocity to distance via a reference wavelength.
source code
float
fnu2flambda(arg, **kwargs)
Switch from Fnu to Flambda via a reference wavelength.
source code
float
flambda2fnu(arg, **kwargs)
Switch from Flambda to Fnu via a reference wavelength.
source code
float
fnu2nufnu(arg, **kwargs)
Switch from Fnu to nuFnu via a reference wavelength.
source code
float
nufnu2fnu(arg, **kwargs)
Switch from nuFnu to Fnu via a reference wavelength.
source code
float
flam2lamflam(arg, **kwargs)
Switch from lamFlam to Flam via a reference wavelength.
source code
float
lamflam2flam(arg, **kwargs)
Switch from lamFlam to Flam via a reference wavelength.
source code
float
distance2spatialfreq(arg, **kwargs)
Switch from distance to spatial frequency via a reference wavelength.
source code
float
spatialfreq2distance(arg, **kwargs)
Switch from spatial frequency to distance via a reference wavelength.
source code
float
per_sr(arg, **kwargs)
Switch from [Q] to [Q]/sr
source code
float
times_sr(arg, **kwargs)
Switch from [Q]/sr to [Q]
source code
float
per_cy(arg, **kwargs)
Switch from radians/s to cycle/s
source code
float
times_cy(arg, **kwargs)
Switch from cycle/s to radians/s
source code
float
period2freq(arg, **kwargs)
Convert period to frequency and back.
source code
 
do_nothing(arg, **kwargs) source code
 
do_sqrt(arg, **kwargs) source code
 
do_quad(arg, **kwargs) source code
 
Hz2_to_ss(change, frequency)
Convert Frequency shift in Hz2 to period change in seconds per second.
source code
    Stellar calibrations
Unit
derive_luminosity(radius, temperature, units=None)
Convert radius and effective temperature to stellar luminosity.
source code
Unit
derive_radius(luminosity, temperature, units='m')
Convert luminosity and effective temperature to stellar radius.
source code
1- or 2-tuple
derive_radius_slo(numax, Deltanu0, teff, unit='Rsol')
Derive stellar radius from solar-like oscillations diagnostics.
source code
 
derive_radius_rot(veq_sini, rot_period, incl=(90.,'deg'), unit='Rsol')
Derive radius from rotation period and inclination angle.
source code
1- or 2-tuple
derive_logg(mass, radius, unit='[cm/s2]')
Convert mass and radius to stellar surface gravity.
source code
1- or 2-tuple
derive_logg_zg(mass, zg, unit='cm s-2', **kwargs)
Convert mass and gravitational redshift to stellar surface gravity.
source code
1- or 2-tuple
derive_logg_slo(teff, numax, unit='[cm/s2]')
Derive stellar surface gravity from solar-like oscillations diagnostics.
source code
 
derive_mass(surface_gravity, radius, unit='kg')
Convert surface gravity and radius to stellar mass.
source code
1- or 2-tuple
derive_zg(mass, logg, unit='cm s-1', **kwargs)
Convert mass and stellar surface gravity to gravitational redshift.
source code
 
derive_numax(mass, radius, temperature, unit='mHz')
Derive the predicted nu_max according to Kjeldsen and Bedding (1995).
source code
 
derive_nmax(mass, radius, temperature)
Derive the predicted n_max according to Kjeldsen and Bedding (1995).
source code
 
derive_Deltanu0(mass, radius, unit='mHz')
Derive the predicted large spacing according to Kjeldsen and Bedding (1995).
source code
 
derive_ampllum(luminosity, mass, temperature, wavelength, unit='ppm')
Derive the luminosity amplitude around nu_max of solar-like oscillations.
source code
 
derive_ampllum_from_velo(velo, temperature, wavelength, unit='ppm')
Derive the luminosity amplitude around nu_max of solar-like oscillations.
source code
 
derive_amplvel(luminosity, mass, unit='cm/s')
Derive the luminosity amplitude around nu_max of solar-like oscillations.
source code
 
derive_galactic_uvw(ra, dec, pmra, pmdec, d, vrad, lsr=False, unit='km s-1')
Calculate the Galactic space velocity (U,V,W) of a star based on the propper motion, location, radial velocity and distance.
source code
 
convert_Z_FeH(Z=None, FeH=None, Zsun=0.0122)
Convert between Z and [Fe/H] using the conversion of Bertelli et al.
source code
    Currencies
 
set_exchange_rates()
Download currency exchange rates from the European Central Bank.
source code
Variables [hide private]
  logger = logging.getLogger("UNITS.CONV")
  _fluxcalib = os.path.join(os.path.abspath(os.path.dirname(__fi...
  _factors = collections.OrderedDict([('m', (1e+00, 'm', 'length...
  _conventions = {'SI': dict(mass= 'kg', length= 'm', time= 's',...
  _names = {'m1': 'length', 's1': 'time', 'kg1': 'mass', 'rad': ...
  _scalings = {'y': 1e-24, 'z': 1e-21, 'a': 1e-18, 'f': 1e-15, '...
  _aliases = [('micron', 'mum'), ('au', 'AU'), ('lbs', 'lb'), ('...
  _switch = {'s1_to_': distance2velocity, 's-1_to_': velocity2di...
Function Details [hide private]

convert(_from, _to, *args, **kwargs)

source code 

Convert one unit to another.

Basic explanation

The unit strings _from and _to should by default be given in the form

erg s-1 cm-2 AA-1

Common alternatives are also accepted (see below).

Square brackets '[]' denote a logarithmic value.

If one positional argument is given, it can be either a scalar, numpy array or uncertainties object. The function will also return one argument of the same type.

If two positional arguments are given, the second argument is assumed to be the uncertainty on the first (scalar or numpy array). The function will also return two arguments.

Basic examples:

>>> convert('km','cm',1.)
100000.0
>>> convert('m/s','km/h',1,0.1)
(3.5999999999999996, 0.36)

Keyword arguments can give extra information, for example when converting from Flambda to Fnu, and should be tuples (float(,error),'unit'):

>>> convert('AA','km/s',4553,0.1,wave=(4552.,0.1,'AA'))
(65.85950307557613, 9.314963362464114)

Extra

The unit strings _from and _to should by default be given in the form

erg s-1 cm-2 AA-1

Common alternatives are also accepted, but don't drive this too far:

erg/s/cm2/AA

The crasiest you're allowed to go is something like

>>> print(convert('10mW m-2/nm','erg s-1 cm-2 AA-1',1.))
1.0

But there is a limit on the interpretation of this prefactor also. Floats will probably not work, and exponentials require exactly two digits.

Parentheses are in no circumstances accepted. Some common aliases are also resolved (for a full list, see dictionary _aliases):

erg/s/cm2/angstrom

You don't really have to spell both units if either the 'from' or 'to' units is consistently within one convention (SI, cgs, solar...). But of course you have to give at least one!:

>>> convert('kg','cgs',1.)
1000.0
>>> convert('g','SI',1.)
0.001
>>> convert('SI','g',1.)
1000.0

WARNINGS:

  1. the conversion involving sr and pixels is not tested.
  2. the conversion involving magnitudes is calibrated but not fully tested
  3. non-integer powers are not functioning yet

Examples:

Spectra:

>>> convert('AA','km/s',4553.,wave=(4552.,'AA'))
65.85950307557613
>>> convert('AA','km/s',4553.,wave=(4552.,0.1,'AA'))
(65.85950307557613, 6.587397133195861)
>>> convert('nm','m/s',455.3,wave=(0.4552,'mum'))
65859.50307564587
>>> convert('km/s','AA',65.859503075576129,wave=(4552.,'AA'))
4553.0
>>> convert('nm','Ghz',1000.)
299792.4579999999
>>> convert('km h-1','nRsol s-1',1.)
0.3993883287866966
>>> print(convert('erg s-1 cm-2 AA-1','SI',1.))
10000000.0

Fluxes:

>>> print(convert('erg/s/cm2/AA','Jy',1e-10,wave=(10000.,'angstrom')))
333.564095198
>>> print(convert('erg/s/cm2/AA','Jy',1e-10,freq=(constants.cc/1e-6,'hz')))
333.564095198
>>> print(convert('erg/s/cm2/AA','Jy',1e-10,freq=(constants.cc,'Mhz')))
333.564095198
>>> print(convert('Jy','erg/s/cm2/AA',333.56409519815202,wave=(10000.,'AA')))
1e-10
>>> print(convert('Jy','erg/s/cm2/AA',333.56409519815202,freq=(constants.cc,'Mhz')))
1e-10
>>> convert('W/m2/mum','erg/s/cm2/AA',1e-10,wave=(10000.,'AA'))
1.0000000000000003e-11
>>> print(convert('Jy','W/m2/Hz',1.))
1e-26
>>> print convert('W/m2/Hz','Jy',1.)
1e+26
>>> print convert('Jy','erg/cm2/s/Hz',1.)
1e-23
>>> print convert('erg/cm2/s/Hz','Jy',1.)
1e+23
>>> print(convert('Jy','erg/s/cm2',1.,wave=(2.,'micron')))
1.49896229e-09
>>> print(convert('erg/s/cm2','Jy',1.,wave=(2.,'micron')))
667128190.396

#>>> print(convert('Jy','erg/s/cm2/micron/sr',1.,wave=(2.,'micron'),ang_diam=(3.,'mas'))) #4511059.82981 #>>> print(convert('Jy','erg/s/cm2/micron/sr',1.,wave=(2.,'micron'),pix=(3.,'mas'))) #3542978.10531 #>>> print(convert('erg/s/cm2/micron/sr','Jy',1.,wave=(2.,'micron'),ang_diam=(3.,'mas'))) #2.21677396826e-07

>>> print(convert('Jy','erg/s/cm2/micron',1.,wave=(2,'micron')))
7.49481145e-10
>>> print(convert('10mW m-2 nm-1','erg s-1 cm-2 AA-1',1.))
1.0

#>>> print convert('Jy','erg s-1 cm-2 micron-1 sr-1',1.,ang_diam=(2.,'mas'),wave=(1.,'micron')) #40599538.4683

Angles:

>>> convert('sr','deg2',1.)
3282.806350011744
>>> ang_diam = 3.21 # mas
>>> scale = convert('mas','sr',ang_diam/2.)
>>> print(ang_diam,scale)
(3.21, 6.054800067947964e-17)
>>> ang_diam = 2*convert('sr','mas',scale)
>>> print(ang_diam,scale)
(3.2100000000000004, 6.054800067947964e-17)

Magnitudes and amplitudes:

>>> print(convert('ABmag','Jy',0.,photband='SDSS.U'))
3767.03798984
>>> print(convert('Jy','erg cm-2 s-1 AA-1',3630.7805477,wave=(1.,'micron')))
1.08848062485e-09
>>> print(convert('ABmag','erg cm-2 s-1 AA-1',0.,wave=(1.,'micron'),photband='SDSS.G'))
4.97510278172e-09
>>> print(convert('erg cm-2 s-1 AA-1','ABmag',1e-8,wave=(1.,'micron'),photband='SDSS.G'))
-0.757994856607
>>> print(convert('ppm','muAmag',1.))
1.0857356618
>>> print(convert('mAmag','ppt',1.,0.1))
(0.9214583192957981, 0.09218827316735488)
>>> print(convert('mag_color','flux_ratio',0.599,0.004,photband='GENEVA.U-B'))
(1.1391327795013377, 0.004196720251233046)

Frequency analysis:

>>> convert('cy/d','muHz',1.)
11.574074074074074
>>> convert('muhz','cy/d',11.574074074074074)
1.0

Interferometry:

>>> convert('m/rad','cy/arcsec',85.,wave=(2.2,'micron'))
187.3143767923207
>>> convert('cm/rad','cy/arcmin',8500.,wave=(2200,'nm'))/60.
187.3143767923207
>>> convert('cy/arcsec','m/rad',187.,wave=(2.2,'mum'))
84.85734129005544
>>> convert('cyc/arcsec','m/rad',187.,wave=(1,'mum'))
38.571518768207014
>>> convert('cycles/arcsec','m',187.,freq=(300000.,'Ghz'))
38.54483473437971
>>> convert('cycles/mas','m',0.187,freq=(300000.,'Ghz'))
38.5448347343797

Temperature:

>>> print(convert('Far','K',123.))
323.705555556
>>> print(convert('kFar','kK',0.123))
0.323705555556
>>> print(convert('K','Far',323.7))
122.99
>>> print(convert('Cel','K',10.))
283.15
>>> print(convert('Cel','Far',10.))
50.0
>>> print(convert('dCel','kFar',100.))
0.05

Time and Dates:

>>> convert('sidereal d','d',1.)
1.0027379093
>>> convert('JD','CD',2446257.81458)
(1985.0, 7.0, 11.314580000005662)
>>> convert('CD','JD',(1985,7,11.31))
2446257.81
>>> convert('CD','JD',(1985,7,11,7,31,59))
2446257.813877315
>>> convert('MJD','CD',0.,jtype='corot')
(2000.0, 1.0, 1.5)
>>> convert('JD','MJD',2400000.5,jtype='mjd')
0.0
>>> convert('MJD','CD',0,jtype='mjd')
(1858.0, 11.0, 17.0)

Coordinates: When converting coordinates with pyephem, you get pyephem coordinates back. They have built-in string represenations and float conversions:

>>> x,y = convert('equatorial','galactic',('17:45:40.4','-29:00:28.1'),epoch='2000')
>>> print (x,y)
(6.282224277178722, -0.000825178833899176)
>>> print x,y
359:56:41.8 -0:02:50.2
>>> x,y = convert('galactic','equatorial',('00:00:00.00','00:00:00.0'),epoch='2000')
>>> print (x,y)
(4.64964430303663, -0.5050315085342665)
>>> print x,y
17:45:37.20 -28:56:10.2

It is also possible to immediately convert to radians or degrees in floats (this can be useful for plotting):

>>> convert('equatorial','deg_coord',('17:45:40.4','-29:00:28.1'))
(266.41833333333335, -29.00780555555556)
>>> convert('equatorial','rad_coord',('17:45:40.4','-29:00:28.1'))
(4.649877104342426, -0.5062817157227474)

Magnetism and Electricity: The stored energy in a magnet, called magnet performance or maximum energy product (often abbreviated BHmax), is typically measured in units of megagauss-oersteds (MGOe). One MGOe is approximately equal to 7957.74715 J/m3 (wikipedia):

>>> convert('MG Oe','J/m3',1.)
7957.747154594768
Parameters:
  • _from (str) - units to convert from
  • _to (str) - units to convert to
  • unpack (boolean, defaults to True) - set to True if you don't want 'uncertainty objects'. If True and uncertainties are given, they will be returned as a tuple (value, error) instead of uncertainty object. Set to False probably only for internal uses
Returns: float
converted value

nconvert(_froms, _tos, *args, **kwargs)

source code 

Convert a list/array/tuple of values with different units to other units.

This silently catches some exceptions and replaces the value with nan!

Returns: array

change_convention(to_, units, origin=None)

source code 

Change units from one convention to another.

Example usage:

>>> units1 = 'kg m2 s-2 K-1 mol-1'
>>> print change_convention('cgs',units1)
K-1 g1 cm2 mol-1 s-2
>>> print change_convention('sol','g cm2 s-2 K-1 mol-1')
Tsol-1 Msol1 Rsol2 mol-1 s-2
Parameters:
  • to_ (str) - convention name to change to
  • units (str) - units to change
  • origin (str) - the name of the original convention
Returns: str
units in new convention

set_convention(units='SI', values='standard', frequency='rad')

source code 

Consistently change the values of the fundamental constants to the paradigm of other system units or programs.

This can be important in pulsation codes, where e.g. the value of the gravitational constant and the value of the solar radius can have siginficant influences on the frequencies.

Setting frequency to rad gives you the following behaviour:

>>> old_settings = set_convention(frequency='rad')
>>> convert('s-1','rad/s',1.0)
1.0
>>> convert('cy/s','rad/s',1.0)
6.283185307179586
>>> convert('Hz','rad/s',1.0)
6.283185307179586
>>> convert('','mas',constants.au/(1000*constants.pc))
0.99999999999623

Setting frequency to cy gives you the following behaviour:

>>> old_settings = set_convention(frequency='cy')
>>> convert('s-1','rad/s',1.0)
6.283185307179586
>>> convert('cy/s','rad/s',1.0)
6.283185307179586
>>> convert('Hz','rad/s',1.0)
6.283185307179586
>>> convert('','mas',constants.au/(1000*constants.pc))
6.2831853071558985
>>> old_settings = set_convention(*old_settings)
Parameters:
  • units (str) - name of the new units base convention
  • values (str) - name of the value set for fundamental constants
Returns: str
name of old convention

get_convention()

source code 

Returns convention and name of set of values.

>>> get_convention()
('SI', 'standard', 'rad')
Returns: str,str
convention, values

get_constant(constant_name, units='SI', value='standard')

source code 

Convenience function to retrieve the value of a constant in a particular system.

>>> Msol = get_constant('Msol','SI')
>>> Msol = get_constant('Msol','kg')
Parameters:
  • constant_name (str) - name of the constant
  • units (str) - name of the unit base system
  • value (str) - name of the parameter set the get the value from
Returns: float
value of the constant in the unit base system

get_constants(units='SI', values='standard')

source code 

Convenience function to retrieve the value of all constants in a particular system.

This can be helpful when you want to attach a copy of the fundamental constants to some module or class instances, and be sure these values never change.

Yes, I know, there's a lot that can go wrong with values that are supposed to be CONSTANT!

Returns: dict

solve_aliases(unit)

source code 

Resolve simple aliases in a unit's name.

Resolves aliases and replaces division signs with negative power.

Parameters:
  • unit (str) - unit (e.g. erg s-1 angstrom-1)
Returns: str
aliases-resolved unit (e.g. erg s-1 AA-1)

components(unit)

source code 

Decompose a unit into: factor, SI base units, power.

You probably want to call solve_aliases first.

Examples:

>>> factor1,units1,power1 = components('m')
>>> factor2,units2,power2 = components('g2')
>>> factor3,units3,power3 = components('hg3')
>>> factor4,units4,power4 = components('Mg4')
>>> factor5,units5,power5 = components('mm')
>>> factor6,units6,power6 = components('W3')
>>> factor7,units7,power7 = components('s-2')
>>> print "{0}, {1}, {2}".format(factor1,units1,power1)
1.0, m, 1
>>> print "{0}, {1}, {2}".format(factor2,units2,power2)
0.001, kg, 2
>>> print "{0}, {1}, {2}".format(factor3,units3,power3)
0.1, kg, 3
>>> print "{0}, {1}, {2}".format(factor4,units4,power4)
1000.0, kg, 4
>>> print "{0}, {1}, {2}".format(factor5,units5,power5)
0.001, m, 1
>>> print "{0}, {1}, {2}".format(factor6,units6,power6)
1.0, m2 kg1 s-3, 3
>>> print "{0}, {1}, {2}".format(factor7,units7,power7)
1.0, s, -2
Parameters:
  • unit (str) - unit name
Returns: (float,str,int)
3-tuple with factor, SI base unit and power

breakdown(unit)

source code 

Decompose a unit into SI base units containing powers.

Examples:

>>> factor1,units1 = breakdown('erg s-1 W2 kg2 cm-2')
>>> factor2,units2 = breakdown('erg s-1 cm-2 AA-1')
>>> factor3,units3 = breakdown('W m-3')
>>> print "{0}, {1}".format(factor1,units1)
0.001, kg5 m4 s-9
>>> print "{0}, {1}".format(factor2,units2)
10000000.0, kg1 m-1 s-3
>>> print "{0}, {1}".format(factor3,units3)
1.0, kg1 m-1 s-3
Parameters:
  • unit (str) - unit's name
Returns: (float,str)
2-tuple factor, unit's base name

compress(unit, ignore_factor=False)

source code 

Compress basic units to more human-readable type.

If you are only interested in a shorter form of the units, regardless of the values, you need to set ignore_factor=True.

For example: erg/s/cm2/AA can be written shorter as W1 m-3, but 1 erg/s/cm2/AA is not equal to 1 W1 m-3. Thus:

>>> print(compress('erg/s/cm2/AA',ignore_factor=True))
W1 m-3
>>> print(compress('erg/s/cm2/AA',ignore_factor=False))
erg/s/cm2/AA

You cannot write the latter shorter without converting the units!

For the same reasoning, compressing non-basic units will not work out:

>>> print(compress('Lsol'))
Lsol
>>> print(compress('Lsol',ignore_factor=True))
W1

So actually, this function is really designed to compress a combination of basic units:

>>> print(compress('kg1 m-1 s-3'))
W1 m-3
Parameters:
  • unit (str) - combination of units to compress
  • ignore_factor (bool) - if True, then a `sloppy compress' will be performed, that is, it is not guaranteed that the multiplication factor equals 1 (see examples)
Returns: str
shorter version of original units

split_units_powers(unit)

source code 

Split unit in a list of units and a list of powers.

>>> y = 'AA A kg-2 F-11 Wb-1'
>>> split_units_powers(y)
(['AA', 'A', 'kg', 'F', 'Wb'], ['1', '1', '-2', '-11', '-1'])
>>> y = 'angstrom A1/kg2 F-11/Wb'
>>> split_units_powers(y)
(['AA', 'A', 'kg', 'F', 'Wb'], ['1', '1', '-2', '-11', '-1'])
>>> y = '10angstrom A1/kg2 F-11/Wb'
>>> split_units_powers(y)
(['10AA', 'A', 'kg', 'F', 'Wb'], ['1', '1', '-2', '-11', '-1'])
>>> y = '10angstrom A1/kg2 5F-11/Wb'
>>> split_units_powers(y)
(['10AA', 'A', 'kg', '5F', 'Wb'], ['1', '1', '-2', '-11', '-1'])
>>> y = '10angstrom A0.5/kg2 5F-11/100Wb2.5'
>>> split_units_powers(y)
(['10AA', 'A', 'kg', '5F', '100Wb'], ['1', '0.5', '-2', '-11', '-2.5'])
Parameters:
  • unit (str) - unit
Returns: list,list
list of unit names, list of powers

is_basic_unit(unit, type)

source code 

Check if a unit is represents one of the basic quantities (length, mass...)

>>> is_basic_unit('K','temperature')
True
>>> is_basic_unit('m','length')
True
>>> is_basic_unit('cm','length')
True
>>> is_basic_unit('km','length') # not a basic unit in any type of convention
False
Parameters:
  • unit (str) - unit to check
  • type (str) - type to check
Returns: bool

is_type(unit, type)

source code 

Check if a unit is of a certain type (mass, length, force...)

>>> is_type('K','temperature')
True
>>> is_type('m','length')
True
>>> is_type('cm','length')
True
>>> is_type('ly','length')
True
>>> is_type('ly','time')
False
>>> is_type('W','power')
True
>>> is_type('kg m2/s3','power')
True
>>> is_type('kg/s2/m','pressure')
True
>>> is_type('W','force')
False
Parameters:
  • unit (str) - unit to check
  • type (str) - type to check
Returns: bool

round_arbitrary(x, base=5)

source code 

Round to an arbitrary base.

Example usage:

>>> round_arbitrary(1.24,0.25)
1.25
>>> round_arbitrary(1.37,0.75)
1.5
Parameters:
  • x (float) - number to round
  • base (integer or float) - base to round to
Returns: float
rounded number

unit2texlabel(unit, full=False)

source code 

Convert a unit string to a nicely formatted TeX label.

For fluxes, also the Flambda, lambda Flambda, or Fnu, or nuFnu is returned.

Parameters:
  • unit (str) - unit name
  • full (bool) - return "name [unit]" or only "unit"
Returns: str

distance2velocity(arg, **kwargs)

source code 

Switch from distance to velocity via a reference wavelength.

Parameters:
  • arg (float) - distance (SI, m)
  • wave (float) - reference wavelength (SI, m)
Returns: float
velocity (SI, m/s)

velocity2distance(arg, **kwargs)

source code 

Switch from velocity to distance via a reference wavelength.

Parameters:
  • arg (float) - velocity (SI, m/s)
  • wave (float) - reference wavelength (SI, m)
Returns: float
distance (SI, m)

fnu2flambda(arg, **kwargs)

source code 

Switch from Fnu to Flambda via a reference wavelength.

Flambda and Fnu are spectral irradiance in wavelength and frequency, respectively

Parameters:
  • arg (float) - spectral irradiance (SI,W/m2/Hz)
  • photband (str ('SYSTEM.FILTER')) - photometric passband
  • wave (float) - reference wavelength (SI, m)
  • freq (float) - reference frequency (SI, Hz)
Returns: float
spectral irradiance (SI, W/m2/m)

flambda2fnu(arg, **kwargs)

source code 

Switch from Flambda to Fnu via a reference wavelength.

Flambda and Fnu are spectral irradiance in wavelength and frequency, respectively

Parameters:
  • arg (float) - spectral irradiance (SI, W/m2/m)
  • photband (str ('SYSTEM.FILTER')) - photometric passband
  • wave (float) - reference wavelength (SI, m)
  • freq (float) - reference frequency (SI, Hz)
Returns: float
spectral irradiance (SI,W/m2/Hz)

fnu2nufnu(arg, **kwargs)

source code 

Switch from Fnu to nuFnu via a reference wavelength.

Flambda and Fnu are spectral irradiance in wavelength and frequency, respectively

Parameters:
  • arg (float) - spectral irradiance (SI,W/m2/Hz)
  • photband (str ('SYSTEM.FILTER')) - photometric passband
  • wave (float) - reference wavelength (SI, m)
  • freq (float) - reference frequency (SI, Hz)
Returns: float
spectral irradiance (SI, W/m2/m)

nufnu2fnu(arg, **kwargs)

source code 

Switch from nuFnu to Fnu via a reference wavelength.

Flambda and Fnu are spectral irradiance in wavelength and frequency, respectively

Parameters:
  • arg (float) - spectral irradiance (SI,W/m2/Hz)
  • photband (str ('SYSTEM.FILTER')) - photometric passband
  • wave (float) - reference wavelength (SI, m)
  • freq (float) - reference frequency (SI, Hz)
Returns: float
spectral irradiance (SI, W/m2/m)

flam2lamflam(arg, **kwargs)

source code 

Switch from lamFlam to Flam via a reference wavelength.

Flambda and Fnu are spectral irradiance in wavelength and frequency, respectively

Parameters:
  • arg (float) - spectral irradiance (SI,W/m2/Hz)
  • photband (str ('SYSTEM.FILTER')) - photometric passband
  • wave (float) - reference wavelength (SI, m)
  • freq (float) - reference frequency (SI, Hz)
Returns: float
spectral irradiance (SI, W/m2/m)

lamflam2flam(arg, **kwargs)

source code 

Switch from lamFlam to Flam via a reference wavelength.

Flambda and Fnu are spectral irradiance in wavelength and frequency, respectively

Parameters:
  • arg (float) - spectral irradiance (SI,W/m2/Hz)
  • photband (str ('SYSTEM.FILTER')) - photometric passband
  • wave (float) - reference wavelength (SI, m)
  • freq (float) - reference frequency (SI, Hz)
Returns: float
spectral irradiance (SI, W/m2/m)

distance2spatialfreq(arg, **kwargs)

source code 

Switch from distance to spatial frequency via a reference wavelength.

Parameters:
  • arg (float) - distance (SI, m)
  • photband (str ('SYSTEM.FILTER')) - photometric passband
  • wave (float) - reference wavelength (SI, m)
  • freq (float) - reference frequency (SI, Hz)
Returns: float
spatial frequency (SI, cy/as)

spatialfreq2distance(arg, **kwargs)

source code 

Switch from spatial frequency to distance via a reference wavelength.

Parameters:
  • arg (float) - spatial frequency (SI, cy/as)
  • photband (str ('SYSTEM.FILTER')) - photometric passband
  • wave (float) - reference wavelength (SI, m)
  • freq (float) - reference frequency (SI, Hz)
Returns: float
distance (SI, m)

per_sr(arg, **kwargs)

source code 

Switch from [Q] to [Q]/sr

Parameters:
  • arg (float) - some SI unit
Returns: float
some SI unit per steradian

times_sr(arg, **kwargs)

source code 

Switch from [Q]/sr to [Q]

Parameters:
  • arg (float) - some SI unit per steradian
Returns: float
some SI unit

Hz2_to_ss(change, frequency)

source code 

Convert Frequency shift in Hz2 to period change in seconds per second.

Frequency in Hz!

derive_luminosity(radius, temperature, units=None)

source code 

Convert radius and effective temperature to stellar luminosity.

Units given to radius and temperature must be understandable by Unit.

Stellar luminosity is returned, by default in Solar units.

I recommend to give the input as follows, since this is most clear also from a programmatical point of view:

>>> print(derive_luminosity((1.,'Rsol'),(5777,'K'),units='erg/s'))
3.83916979644e+33 erg/s

The function is pretty smart, however: it knows what the output units should be, so you could ask to put them in the 'standard' units of some convention:

>>> print(derive_luminosity((1.,'Rsol'),(5777,'K'),units='cgs'))
3.83916979644e+33 g1 cm2 s-3

If you give nothing, everything will be interpreted in the current convention's units:

>>> print(derive_luminosity(7e8,5777.))
3.88892117726e+26 W1

You are also free to give units (you can do this in a number of ways), and if you index the return value, you can get the value ([0], float or uncertainty) or unit ([1],string):

>>> derive_luminosity((1.,'Rsol'),(1.,0.1,'Tsol'),units='Lsol')[0]
1.0000048246799305+/-0.4000019298719723
>>> derive_luminosity((1.,'Rsol'),((1.,0.1),'Tsol'),units='Lsol')[0]
1.0000048246799305+/-0.4000019298719723

Finally, you can simply work with Units:

>>> radius = Unit(1.,'Rsol')
>>> temperature = Unit(5777.,'K')
>>> print(derive_luminosity(radius,temperature,units='foe/s'))
3.83916979644e-18 foe/s

Everything should be independent of the current convention, unless it needs to be:

>>> old = set_convention('cgs')
>>> derive_luminosity((1.,'Rsol'),(1.,0.1,'Tsol'),units='Lsol')[0]
1.0000048246799305+/-0.4000019298719722
>>> print(derive_luminosity(7e10,5777.))
3.88892117726e+33 erg1 s-1
>>> sol = set_convention('sol')
>>> print(derive_luminosity(1.,1.))
3.9982620567e-22 Msol1 Rsol2 s-3
>>> print(derive_luminosity(1.,1.,units='Lsol'))
1.00000482468 Lsol
>>> old = set_convention(*old)
Parameters:
  • radius (2 or 3 tuple, Unit or float (current convention)) - (radius(, error), units)
  • temperature (2 or 3 tuple, Unit or float (current convention)) - (effective temperature(, error), units)
Returns: Unit
luminosity

derive_radius(luminosity, temperature, units='m')

source code 

Convert luminosity and effective temperature to stellar radius.

Units given to luminosity and temperature must be understandable by convert.

Stellar radius is returned in SI units.

Example usage:

>>> print(derive_radius((3.9,'[Lsol]'),(3.72,'[K]'),units='Rsol'))
108.091293736 Rsol
>>> print(derive_radius(3.8e26,5777.,units='Rjup'))
9.89759670363 Rjup
Parameters:
  • luminosity (2 or 3 tuple) - (Luminosity(, error), units)
  • temperature (2 or 3 tuple, Unit or float (current convention)) - (effective temperature(, error), units)
Returns: Unit
radius

derive_radius_slo(numax, Deltanu0, teff, unit='Rsol')

source code 

Derive stellar radius from solar-like oscillations diagnostics.

Large separation is frequency spacing between modes of different radial order.

Small separation is spacing between modes of even and odd angular degree but same radial order.

Parameters:
  • numax (2 or 3 tuple) - (numax(, error), units)
  • Deltanu0 (2 or 3 tuple) - (large separation(, error), units)
  • teff (2 or 3 tuple) - (effective temperature(, error), units)
Returns: 1- or 2-tuple
radius (and error) in whatever units

derive_logg(mass, radius, unit='[cm/s2]')

source code 

Convert mass and radius to stellar surface gravity.

Units given to mass and radius must be understandable by convert.

Logarithm of surface gravity is returned in CGS units.

Parameters:
  • mass (2 or 3 tuple) - (mass(, error), units)
  • radius (2 or 3 tuple) - (radius(, error), units)
Returns: 1- or 2-tuple
log g (and error) in CGS units

derive_logg_zg(mass, zg, unit='cm s-2', **kwargs)

source code 

Convert mass and gravitational redshift to stellar surface gravity. Provide the gravitational redshift as a velocity.

Units given to mass and zg must be understandable by convert.

Logarithm of surface gravity is returned in CGS units unless otherwhise stated in unit.

>>> print derive_logg_zg((0.47, 'Msol'), (2.57, 'km/s'))
5.97849814386
Parameters:
  • mass (2 or 3 tuple) - (mass(, error), units)
  • zg (2 or 3 tuple) - (zg(, error), units)
  • unit (str) - unit of logg
Returns: 1- or 2-tuple
log g (and error)

derive_logg_slo(teff, numax, unit='[cm/s2]')

source code 

Derive stellar surface gravity from solar-like oscillations diagnostics.

Units given to teff and numax must be understandable by convert.

Logarithm of surface gravity is returned in CGS units.

Parameters:
  • teff (2 or 3 tuple) - (effective temperature(, error), units)
  • numax (2 or 3 tuple) - (numax(, error), units)
Returns: 1- or 2-tuple
log g (and error) in CGS units

derive_zg(mass, logg, unit='cm s-1', **kwargs)

source code 

Convert mass and stellar surface gravity to gravitational redshift.

Units given to mass and logg must be understandable by convert.

Gravitational redshift is returned in CGS units unless otherwhise stated in unit.

>>> print derive_zg((0.47, 'Msol'), (5.98, 'cm s-2'), unit='km s-1')
2.57444756874
Parameters:
  • mass (2 or 3 tuple) - (mass(, error), units)
  • logg (2 or 3 tuple) - (logg(, error), units)
  • unit (str) - unit of zg
Returns: 1- or 2-tuple
log g (and error)

derive_numax(mass, radius, temperature, unit='mHz')

source code 

Derive the predicted nu_max according to Kjeldsen and Bedding (1995).

Example: compute the predicted numax for the Sun in mHz >>> print derive_numax((1.,'Msol'),(1.,'Rsol'),(5777.,'K')) 3.05

derive_nmax(mass, radius, temperature)

source code 

Derive the predicted n_max according to Kjeldsen and Bedding (1995).

Example: compute the predicted numax for the Sun in mHz >>> print derive_nmax((1.,'Msol'),(1.,'Rsol'),(5777.,'K')) 21.0

derive_Deltanu0(mass, radius, unit='mHz')

source code 

Derive the predicted large spacing according to Kjeldsen and Bedding (1995).

Example: compute the predicted large spacing for the Sun in mHz >>> print derive_Deltanu0((1.,'Msol'),(1.,'Rsol'),unit='muHz') 134.9

derive_ampllum(luminosity, mass, temperature, wavelength, unit='ppm')

source code 

Derive the luminosity amplitude around nu_max of solar-like oscillations.

See Kjeldsen and Bedding (1995).

>>> print derive_ampllum((1.,'Lsol'),(1,'Msol'),(5777.,'K'),(550,'nm'))
(4.7, 0.3)

derive_ampllum_from_velo(velo, temperature, wavelength, unit='ppm')

source code 

Derive the luminosity amplitude around nu_max of solar-like oscillations.

See Kjeldsen and Bedding (1995).

>>> print derive_ampllum_from_velo((1.,'m/s'),(5777.,'K'),(550,'nm'))
(20.1, 0.05)

derive_amplvel(luminosity, mass, unit='cm/s')

source code 

Derive the luminosity amplitude around nu_max of solar-like oscillations.

See Kjeldsen and Bedding (1995).

>>> print derive_amplvel((1.,'Lsol'),(1,'Msol'),unit='cm/s')
(23.4, 1.4)

derive_galactic_uvw(ra, dec, pmra, pmdec, d, vrad, lsr=False, unit='km s-1')

source code 

Calculate the Galactic space velocity (U,V,W) of a star based on the propper motion, location, radial velocity and distance. (U,V,W) are returned in km/s unless stated otherwhise in unit.

Follows the general outline of Johnson & Soderblom (1987, AJ, 93,864) except that U is positive outward toward the Galactic *anti*center, and the J2000 transformation matrix to Galactic coordinates is taken from the introduction to the Hipparcos catalog.

Uses solar motion from Coskunoglu et al. 2011 MNRAS: (U,V,W)_Sun = (-8.5, 13.38, 6.49)

Example for HD 6755:

>>> derive_galactic_uvw((17.42943586, 'deg'), (61.54727506, 'deg'), (628.42, 'mas yr-1'), (76.65, 'mas yr-1'), (139, 'pc'), (-321.4, 'km s-1'), lsr=True)
(142.66328352779027, -483.55149105148121, 93.216106970932813)
Parameters:
  • ra - Right assension of the star
  • dec - Declination of the star
  • pmra - propper motion in RA
  • pmdec - propper motion in DEC
  • d - distance
  • vrad - radial velocity
  • lsr - if True, correct for solar motion
  • unit - units for U,V and W
Returns:
(U,V,W)

convert_Z_FeH(Z=None, FeH=None, Zsun=0.0122)

source code 

Convert between Z and [Fe/H] using the conversion of Bertelli et al. (1994). Provide one of the 2 arguments, and the other will be returned.

log10(Z/Zsun) = 0.977 [Fe/H]

set_exchange_rates()

source code 

Download currency exchange rates from the European Central Bank.

Decorators:
  • @memoized

Variables Details [hide private]

_fluxcalib

Value:
os.path.join(os.path.abspath(os.path.dirname(__file__)), 'fluxcalib.da\
t')

_factors

Value:
collections.OrderedDict([('m', (1e+00, 'm', 'length', 'meter')), ('AA'\
, (1e-10, 'm', 'length', 'angstrom')), ('AU', (constants.au, constants\
.au_units, 'length', 'astronomical unit')), ('pc', (constants.pc, cons\
tants.pc_units, 'length', 'parsec')), ('ly', (constants.ly, constants.\
ly_units, 'length', 'light year')), ('bs', (5e-9, 'm', 'length', 'bear\
d second')), ('Rsol', (constants.Rsol, constants.Rsol_units, 'length',\
 'Solar radius')), ('Rearth', (constants.Rearth, constants.Rearth_unit\
s, 'length', 'Earth radius')), ('Rjup', (constants.Rjup, constants.Rju\
...

_conventions

Value:
{'SI': dict(mass= 'kg', length= 'm', time= 's', temperature= 'K', elec\
tric_current= 'A', lum_intens= 'cd', amount= 'mol', currency= 'EUR'), \
'cgs': dict(mass= 'g', length= 'cm', time= 's', temperature= 'K', elec\
tric_current= 'A', lum_intens= 'cd', amount= 'mol', currency= 'EUR'), \
'sol': dict(mass= 'Msol', length= 'Rsol', time= 's', temperature= 'Tso\
l', electric_current= 'A', lum_intens= 'cd', amount= 'mol', currency= \
'EUR'), 'imperial': dict(mass= 'lb', length= 'yd', time= 's', temperat\
ure= 'K', electric_current= 'A', lum_intens= 'cd', amount= 'mol', curr\
...

_names

Value:
{'m1': 'length', 's1': 'time', 'kg1': 'mass', 'rad': 'angle', 'deg': '\
angle', 'K': 'temperature', 'm1 s-2': 'acceleration', 'kg1 s-2': 'surf\
ace tension', 'cy-1 kg1 s-2': 'flux density', 'kg1 m-1 s-3': 'flux den\
sity', 'kg1 s-3': 'flux', 'cy1 s-1': 'frequency', 'cy-1 kg1 rad-2 s-2'\
: 'specific intensity', 'kg1 m-1 rad-2 s-3': 'specific intensity', 'kg\
1 rad-2 s-3': 'total intensity', 'kg1 m2 s-3': 'luminosity', 'm1 s-1':\
 'velocity', 'kg1 m-1 s-2': 'pressure', 'm2': 'area', 'm3': 'volume', \
'kg1 m2 s-2': 'energy', 'kg1 m1 s-2': 'force', 'A': 'electric current'\
...

_scalings

Value:
{'y': 1e-24, 'z': 1e-21, 'a': 1e-18, 'f': 1e-15, 'p': 1e-12, 'n': 1e-0\
9, 'mu': 1e-06, 'u': 1e-06, 'm': 1e-03, 'c': 1e-02, 'd': 1e-01, 'da': \
1e+01, 'h': 1e+02, 'k': 1e+03, 'M': 1e+06, 'G': 1e+09, 'T': 1e+12, 'P'\
: 1e+15, 'E': 1e+18, 'Z': 1e+21, 'Y': 1e+24}

_aliases

Value:
[('micron', 'mum'), ('au', 'AU'), ('lbs', 'lb'), ('micro', 'mu'), ('mi\
lli', 'm'), ('kilo', 'k'), ('mega', 'M'), ('giga', 'G'), ('nano', 'n')\
, ('dyne', 'dyn'), ('watt', 'W'), ('Watt', 'W'), ('Hz', 'hz'), ('liter\
', 'l'), ('litre', 'l'), ('Tesla', 'T'), ('Ampere', 'A'), ('Coulomb', \
'C'), ('joule', 'J'), ('Joule', 'J'), ('jansky', 'Jy'), ('Jansky', 'Jy\
'), ('jy', 'Jy'), ('arcsec', 'as'), ('arcmin', 'am'), ('cycles', 'cy')\
, ('cycle', 'cy'), ('cyc', 'cy'), ('angstrom', 'AA'), ('Angstrom', 'AA\
'), ('inch', 'in'), ('stone', 'st'), ('^', ''), ('**', ''), ('celcius'\
...

_switch

Value:
{'s1_to_': distance2velocity, 's-1_to_': velocity2distance, 'm1rad-1_t\
o_': distance2spatialfreq, 'm-1rad1_to_': spatialfreq2distance, 'm1rad\
1s1_to_': fnu2flambda, 'm-1rad-1s-1_to_': flambda2fnu, 'm1rad-1s1_to_'\
: fnu2flambda, 'm-1rad1s-1_to_': flambda2fnu, 'rad-1s1_to_': fnu2nufnu\
, 'rad1s-1_to_': nufnu2fnu, 'm1_to_': lamflam2flam, 'm-1_to_': flam2la\
mflam, 'sr1_to_': per_sr, 'sr-1_to_': times_sr, 'rad1_to_': do_nothing\
, 'rad-1_to_': do_nothing, 'rad-2sr1_to_': do_nothing, 'rad2sr-1_to_':\
 do_nothing, 'rad-1sr1_to_': do_sqrt, 'rad1sr-1_to_': do_quad, 'rad-1s\
...