Magnitude conversions

Translate absolute magnitude to power scale

PyAstronomy.pyasl.absMagToPower(am, absMagSun=4.75, absLumSun=3.846e+33)

Convert absolute magnitude to power scale

The default values for the absolute magnitude and luminosity of the Sun are adopted from Harmanec and Prsa 2011 (2011PASP..123..976H).

Parameters
amfloat

Absolute magnitude.

absMagSunfloat, optional

Absolute magnitude of the Sun.

absLumSunfloat, optional

Absolute luminosity of the Sun. The default is given in units of erg/s.

Returns
Powerfloat

Total emitted power. Same units as absLumSun; the default corresponds to erg/s.

Example

from __future__ import print_function, division
from PyAstronomy import pyasl

absMagSun = 4.75
print("Absolute bolometric magnitude of the Sun: ", absMagSun)
print("  Absolute luminosity [erg/s]: ", pyasl.absMagToPower(absMagSun))

Translate distance module into distance

PyAstronomy.pyasl.absModuleToDist(magApp, magAbs)

Convert apparent and absolute magnitude into distance.

Parameters
magAppfloat

Apparent magnitude of object.

magAbsfloat

Absolute magnitude of object.

Returns
Distancefloat

The distance resulting from the difference in apparent and absolute magnitude [pc].

Example

from __future__ import print_function, division
from PyAstronomy import pyasl

# Apparent magnitude
appMag = 11.37
# Absolute (bolometric) magnitude of Sun
absMagSun = 4.75

print("Distance of a sun-like star with apparent bolometric ", end=' ')
print("brightness of 11.37 mag: %5.2f pc" % (pyasl.absModuleToDist(appMag, absMagSun)))

Convert magnitude into flux density

PyAstronomy.pyasl.magToFluxDensity_bessel98(band, mag, mode='nu')

Convert magnitude into flux density according to Bessel et al. 1998

The conversion implemented here is based on the data given in Table A2 of Bessel et al. 1998, A&A 333, 231-250, which gives “Effective wavelengths (for an A0 star), absolute fluxes (corresponding to zero magnitude) and zeropoint magnitudes for the UBVRI- JHKL Cousins-Glass-Johnson system”. Note that zp(f_nu) and zp(f_lam) are exchanged in the original table.

Parameters
bandstring

Any of U, B, V, R, I, J, H, K, Kp, L, and L*

magfloat, array

The magnitude value to be converted

modestring, {nu, mod}

Determines whether f_nu or f_lam will be calculated.

Returns
f_nu/lamfloat

The corresponding flux density in units if erg/cm**2/s/Hz in the case of mode ‘nu’ and erg/cm**2/s/A in the case of ‘lam’.

lam_efffloat

Effective filter wavelength in Angstrom

Example

from __future__ import print_function
from PyAstronomy import pyasl
import numpy as np

mag_R = 15.5

fd_nu, le = pyasl.magToFluxDensity_bessel98("R", mag_R, "nu")
fd_lam, _ = pyasl.magToFluxDensity_bessel98("R", mag_R, "lam")

print("R-band magnitude: ", mag_R)
print("R-band flux density [erg/cm**2/s/Hz]: ", fd_nu)
print("R-band flux density [erg/cm**2/s/A]: ", fd_lam)

print("Effective wavelength of filter [A]: ", le)
print("Convert f_nu into f_lam [erg/cm**2/s/A] by multiplication with (c/lam**2): ",
      fd_nu * (299792458e2/(le/1e8)**2) / 1e8)