irfpy.util.solarspectra

Solar spectrum models at the Earth orbit

Code author: Yoshifumi Futaana

BlackBody([temperature])

Black body assumption of the Sun.

BlackBodyModFuv(*args, **kwds)

Black body assumption with constant FUV (<~100 nm) modification

AM0()

Air mass 0 model (at extraterrestorial) under ASTM E-490.

Solar spectrum models at the Earth orbit (1 AU). See also irfpy.util.planck.

class irfpy.util.solarspectra.BlackBody(temperature=5778.0)[source]

Bases: object

Black body assumption of the Sun.

Sun is proxied by a black body with tempereature 5778 K.

Using the planck’s formula (see irfpy.util.planck.BlackBodyMksa) the spectrum at the Sun, i.e. power by area by wavelength [W / m2 / sr/ m], is immediately known as \(B(\lambda, T)\). This is implemented in irfpy.util.planck.

The total spectral power is the integral over the Sun surface area and the emission angle (in solid angle). As isotropic emission is assumed, the solid angle integral can be done

\[\int{\cos\theta d\Omega} = \int_0^{2\pi}d\phi\int_0^{\pi/2}\cos\theta\sin\theta d\theta = \pi\]

The \(\cos\theta\) appears from the inner product of the flux and the surface area. See also how to derive Stefan-Boltzmann law. Same calculation is done.

Anyway, the total spectral power, in the unit of [W / m], is

\[F(\lambda, T) = \pi B(\lambda, T) \times 4\pi Rs^2\]

Here Rs is the radius of the Sun.

Consider a unit area at the Earth’s orbit. The spectral power through the area is, as the emission is isotropic, calculated by multiplying the fraction of the area to the sphare surface at 1AU, namely,

\[u(\lambda, T) = F(\lambda, T) \times \frac{1}{4\pi a^2}\]

where a is 1AU. The unit of u is [W / m2 / m].

In summary, the spectrum power at the Earth orbit, u [W / m2 / m] is

\[u(\lambda, T) = \pi \frac{Rs^2}{a^2} B(\lambda, T)\]

Substituting Rs=695500 km and a=149597871 km,

\[u(\lambda, T) = 6.8687\times10^5 B(\lambda, T)\]
name

Name of the model

spectral_irradiance(wavelength)[source]

Return the spectral irradiance (W / m2 / m).

Return the spectral irradiance, in the unit of W / m2 / m.

Parameters:

wavelength – Wavelength in m.

Returns:

Spectral irradiance in W / m2 / m

>>> sun = BlackBody(5800)    # 5800 K blackbody approx.
>>> irr500nm = sun.spectral_irradiance(500e-9)
>>> print('%.3e' % irr500nm)
1.846e+09
photon_flux(wavelength)[source]

Return the photon’s differential flux (# / m2 / s / m).

The spectral irradiance, u [W / m2 / m = J / m2 / s / m], can be converted to the photon flux.

A photon with wavelength, \(\lambda\), will have the energy \(E = hc/\lambda\). Thus, the number of photons through a unit area per unit time is \(F = u / E\).

>>> sun = BlackBody(5800)
>>> fph = sun.photon_flux(500e-9)   # in /m2 /s /m
>>> fph = fph / 1e9   # in m2 /s /nm
>>> print('%.3e' % fph)
4.648e+18
class irfpy.util.solarspectra.BlackBodyModFuv(*args, **kwds)[source]

Bases: BlackBody

Black body assumption with constant FUV (<~100 nm) modification

name

Name of the model

spectral_irradiance(wavelength)[source]
>>> sun = BlackBodyModFuv(5800)
>>> print('%.3e' % (sun.spectral_irradiance(500e-9)))
1.846e+09
class irfpy.util.solarspectra.AM0[source]

Bases: object

Air mass 0 model (at extraterrestorial) under ASTM E-490.

>>> am0 = AM0()
>>> irr = am0.spectral_irradiance(500e-9)
>>> print('%.3e' % irr)
1.914e+09
>>> fph = am0.photon_flux(np.array([500e-9, 1000e-9]))
>>> print('%.3e %.3e' % (fph[0] / 1e9, fph[1] / 1e9))   # in m2 /nm /s
4.816e+18 3.765e+18
logger = <Logger irfpy.util.solarspectra.AM0 (DEBUG)>
irradfunc

The irradiance function from wavelength (nm) to irradiance (W/m2/nm)

spectral_irradiance(wavelength)[source]

Return the irradiance with interface of MKSA

photon_flux(wavelength)[source]
irfpy.util.solarspectra.doctests()[source]