Functions for blackbody thermal emission (and its derivatives) calculation

PlanckFunctions.PlanckFunctionsModule

PlanckFunctions module provides a set of functions for evaluating the Planck thermal emission spectrum intensity (spectral radiance) and its derivatives with respect to wavelength and temperature. It also provides function to evaluate the integral over the wavelength radiance and Rosseland-averaged and Planck-averaged of spectral coefficients.

Wavelength units are microns and all temperatures should be in Kelvins.

Main functions are: ibb - spectral intensity (spectral radiance) in [W/m²⋅sr⋅μm]

∇ₜibb - spectral intensity first derivative with respect to temperature

∇²ₜibb - spectral intensity second derivative with respect to temperature

∇ₗibb - spectral intensity first derivative with respect to wavelength

∇²ₗibb - spectral intensity second derivative with respect to wavelength

rosseland_averaged_attenuation - Rosseland-averaging of spectral attenuation

planck_averaged_attenuation - Planck-averaging of spectral attenuation

Main literature sources are:

J.R.Howell,M.P.Menguc,J.R.Howell,M.P.Menguc,K.Daun,R.Siegel. Thermal radiation heat transfer. Seventh edition. 2021

Risch, T.K., User's Manual: Routines for Radiative Heat Transfer and Thermometry. NASA/TM-2016-219103. 2016, Edwards, California: Armstrong flight Research Center

source
PlanckFunctions.C₁Constant
C₁ constant for Planck function multiplier in [W⋅μm⁴/(m²⋅sr)]
source J.R.Howell,M.P.Menguc,J.R.Howell,M.P.Menguc,K.Daun,R.Siegel. Thermal radiation heat transfer. Seventh edition. 2021  see Appendix A.
Multiplied by 2 with respect to the source.
source
PlanckFunctions.C₂Constant
C₂ constant for Planck spectral intensity exponent in [μm⋅K]
source J.R.Howell,M.P.Menguc,J.R.Howell,M.P.Menguc,K.Daun,R.Siegel. Thermal radiation heat transfer. Seventh edition. 2021  see Appendix A.
source
PlanckFunctions.C₃Constant
C₃ constant of Wien's displacement law [μm⋅K]
source J.R.Howell,M.P.Menguc,J.R.Howell,M.P.Menguc,K.Daun,R.Siegel. Thermal radiation heat transfer. Seventh edition. 2021  see Appendix A.
source
PlanckFunctions.C₄Constant
C₄ constant in equation for maximum blackbody intensity [W/(m²⋅μm⋅sr*K⁵)]
source J.R.Howell,M.P.Menguc,J.R.Howell,M.P.Menguc,K.Daun,R.Siegel. Thermal radiation heat transfer. Seventh edition. 2021  see Appendix A.
source
PlanckFunctions.σConstant
Stefan-Boltzmann constant [W/(m²*K⁴)]
source J.R.Howell,M.P.Menguc,J.R.Howell,M.P.Menguc,K.Daun,R.Siegel. Thermal radiation heat transfer. Seventh edition. 2021  see Appendix A.
source
PlanckFunctions.DₗibbMethod
Dₗibb(λ,T)

Returns a three-element tuple of (1.bb intensity,2.its first and 3.second derivative 
with respect to the wavelentgh)

Input:
    λ - wavelength, μm
    T - temperature, K
source
PlanckFunctions.Dₜibb!Method
Dₜibb!(input_tuple, λ::AbstractVector,T)

In-place filling the tuple of (bb intensity, its first ,and second ) derivatives with 
respect to temperature
d²Ibb/dT² = C₁*(eᵃ¹/(eᵃ¹-1)²)*(C₂/(λ⁶*T³))*[(C₂/(λ*T))*(2*eᵃ¹/(eᵃ¹-1)-1)-2]
d²Ibb/dT² = C₁*(eᵃ¹/(eᵃ¹-1)²)*(a₁/(λ⁵*T²))*[a₁*(2*eᵃ¹/(eᵃ¹-1) -1)-2]
a₁=C₂/(λ*T)
a₂ = 1/(eᵃ¹-1)   #  1/expm1(a1)
a₃ = eᵃ¹/(eᵃ¹-1) #  exp(a)/expm1(a)
d²Ibb/dT² = C₁*a₂*a₃*(a₁/(λ⁵*T²))*[a₁*(2*a₃ - 1))-2]
    as far as 
        Ibb = (λ⁻⁵)* C₁*a₂
    and 
        dIbb/dT = C₁*a₃*a₂*a₁*(1/(λ⁵*T)) = a₃*a₁*Ibb/T 
    hense
        d²Ibb/dT² = C₁*a₂*a₃*a₁*(1/(λ⁵*T²))*[a₁*(2*a₃ - 1))-2] 
            = [a₃*a₁*Ibb/T^2]*[a₁*(2*a₃ - 1))-2] 
                = [(dIbb/dT)/T]*[a₁*(2*a₃ - 1))-2] 
Input:
    input_tuple, [Nx0 vector or nothing,Nx0 vector or nothing, Nx0 vector or nothing]
    λ - wavelength, μm, [Nx0]
    T - temperature, K
source
PlanckFunctions.DₜibbMethod
Dₜibb(λ::AbstractVector,T::AbstractVector)

Calculates tuple of (Ibb,dIbb/dT,d²Ibb/dT²) calculated according to:
d²Ibb/dT² = C₁*(eᵃ¹/(eᵃ¹-1)²)*(C₂/(λ⁶*T³))*[(C₂/(λ*T))*(2*eᵃ¹/(eᵃ¹-1)-1)-2]
d²Ibb/dT² = C₁*(eᵃ¹/(eᵃ¹-1)²)*(a₁/(λ⁵*T²))*[a₁*(2*eᵃ¹/(eᵃ¹-1) -1)-2]
a₁=C₂/(λ*T)
a₂ = 1/(eᵃ¹-1)   #  1/expm1(a1)
a₃ = eᵃ¹/(eᵃ¹-1) #  exp(a)/expm1(a)
d²Ibb/dT² = C₁*a₂*a₃*(a₁/(λ⁵*T²))*[a₁*(2*a₃ - 1))-2]
    as far as 
        Ibb = (λ⁻⁵)* C₁*a₂
    and 
        dIbb/dT = C₁*a₃*a₂*a₁*(1/(λ⁵*T)) = a₃*a₁*Ibb/T 
    hense
        d²Ibb/dT² = C₁*a₂*a₃*a₁*(1/(λ⁵*T²))*[a₁*(2*a₃ - 1))-2] 
            = [a₃*a₁*Ibb/T^2]*[a₁*(2*a₃ - 1))-2] 
                = [(dIbb/dT)/T]*[a₁*(2*a₃ - 1))-2] 
Input:
    λ - wavelength region, μm
    T - temperature, Kelvins
Returns:
    (Ibb,dIbb/dT,d²Ibb/dT²)
source
PlanckFunctions.a₁₂₃!Method
a₁₂₃!(amat::AbstractMatrix,λ::AbstractVector,T::Float64)

In-place filling of the intermediate matrix
a₁=C₂/(λ*T)  - amat first column
a₂ = 1/(eᵃ¹-1)  - amat second column 
a₃ = eᵃ¹/(eᵃ¹-1) - amat third column

Input:
    amat - matrix of intermediate coefficients size [Nx3]
    λ - wavelength in μm,  [Nx0]
    T - temperature in Kelvins
source
PlanckFunctions.a₁₂₃Method
a₁₂₃(λ::Float64,T::Float64)

Return tuple of three values:

1 - a1 = C₂/(λ*T)
2 - a2 = 1/expm1(a1)
3 - a3 = exp(a1)*a2
source
PlanckFunctions.band_powerMethod
band_power(T;λₗ=0.0,λᵣ=Inf,tol=1e-6)

Total bb with temperature T integral intensity within (in-band radiance), [W/(m²⋅sr)]
the spectral range λₗ...λᵣ (by default the range is 0...inf)
tol - tolerance of intehration

Input:
    T - temperature,Kelvins
    (optional)
    λₗ - left wavelength boundary, μm
    λᵣ - right wavelength boundary, μm
    tol - intergation tolerance
source
PlanckFunctions.ibb!Method
ibb!(i::AbstractVector,λ::AbstractVector,amat::AbstractMatrix)::Nothing

In-place blackbody intensity with intermediate coefficients provided externally, [W/m2-sr-mkm]
Ibb =  C₁*(λ⁻⁵)*a₂ , where
a₁=C₂/(λ*T)  - amat first column
a₂ = 1/(eᵃ¹-1)  - amat second column 

Input:
    i - BB intensity, [Nx0]
    λ - wavelength in μm,  [Nx0]
    amat - matrix of intermediate coefficients,  [Nx3]
source
PlanckFunctions.ibb!Method
ibb!(i::AbstractVector,λ::AbstractVector,T::Float64)

In-place blackbody intensity,  [W/m²⋅sr⋅μm]
Ibb = (λ⁻⁵)* C₁/(eᵃ¹-1) , where a₁=C₂/(λ*T)
Input:
    i - bb intensity vector, [Nx0]
    λ - wavelength in μm, [Nx0]
    T - temperature in Kelvins
source
PlanckFunctions.ibbMethod
ibb(λ::AbstractVector,amat::AbstractMatrix)

Blackbody spectral intensity (spectral radiance)  with intermediate matrix provided externally, [W/m²⋅sr⋅μm]
Ibb =  C₁*(λ⁻⁵)*a₂ , where
a₁=C₂/(λ*T)  - amat first column
a₂ = 1/(eᵃ¹-1)  - amat second column 
Input:
    amat - matrix of intermediate coefficients,  [Nx3]
    λ - wavelength in μm,  [Nx0]
source
PlanckFunctions.ibbMethod
ibb(λ::AbstractVector,T::AbstractVector)

Blackbody spectral intensity (spectral radiance),  [W/m²⋅sr⋅μm]
Ibb = (λ⁻⁵)* C₁/(eᵃ¹-1) , where a₁=C₂/(λ*T)
Input:
    λ - wavelength in μm, [Nx0]
    T - temperature in Kelvins [Mx0]
source
PlanckFunctions.ibbMethod
ibb(λ,T)

Blackbody spectral intensity (spectral radiance), [W/m²⋅sr⋅μm]
Ibb = (λ⁻⁵)* C₁/(eᵃ¹-1) , where a₁=C₂/(λ*T)
Input:
    λ - wavelength in μm
    T - temperature in Kelvins
source
PlanckFunctions.planck_averagedMethod
planck_averaged(x::AbstractVector, λ::AbstractVector,T::Number)

Evaluates the Planck-averaged value of x(λ) for temperature T:

xᵣ = ∫x(λ)ibb(λ,T)dλ/∫ibb(λ,T)dλ 

Can be used for example to evaluate the integral emissiovity from measured spectral emissivity
source
PlanckFunctions.planck_averaged_attenuationMethod
planck_averaged_attenuation(α::AbstractVector, λ::AbstractVector,T::Number)

Evaluates the Planck-averaged spectral attenuation coefficient (the summation of
spectral scattering and absorption coefficients) α(λ) for temperature T:

αᵣ = (∫(1/α(λ))⋅ibb(λ,T)dλ/∫ibb(λ,T)dλ)⁻¹
source
PlanckFunctions.powerMethod
power(T)

Returns integral (over the wavelength) intensity of BB (radiance)  at temperature T

Units: W/(m²⋅sr)

Input:
    T - temperature, K
source
PlanckFunctions.rosseland_averaged_attenuationMethod
rosseland_averaged_attenuation(α::AbstractVector, λ::AbstractVector,T::Number)

Evaluates the Rosseland-averaged spectral attenuation coefficient (the summation of
spectral scattering and absorption coefficients) α(λ) for temperature T:

αᵣ = (∫(1/α(λ))⋅∇ₜibb(λ,T)dλ/∫∇ₜibb(λ,T)dλ)⁻¹
source
PlanckFunctions.weighted_averageFunction
weighted_average(α::AbstractVector, λ::AbstractVector,T,g::Union{typeof(ibb),typeof(∇ₜibb),typeof(∇²ₜibb)},f::Function=identity)

Generic function to evaluate the averaged value of some `f(x)` function of variable `x` dependent
on wavelength `λ` for temperature `T`. Uses linear approximation for the discrete variable and square polynomial 
for the g function 

xᵣ = ∫f(x)g(λ,T)dλ/∫g(λ,T)dλ, the default value of f is inv, e.g. if f = inv:
xᵣ = ∫g(λ,T)/x(λ)dλ/∫g(λ,T)dλ
source
PlanckFunctions.λₘMethod
λₘ(T)

The wavelength (in μm) of bb intensity maximum vs temperature T 
argmax(Planck(T))  - Wien's displacement law

Input:
    T - temperature in Kelvins
source
PlanckFunctions.∇²ₗibbMethod
∇²ₗibb(λ,T)

BB intensity second derivative with respect to the wavelength

Input:
    λ - wavelength, μm
    T - temperature, K
source
PlanckFunctions.∇²ₜibb!Method
∇²ₜibb!(h::AbstractMatrix{Float64},λ::AbstractVector{Float64},T::AbstractVector{Float64})

In-place bb intensity second order derivative with respect to temperature

d²Ibb/dT² = C₁*(eᵃ¹/(eᵃ¹-1)²)*(C₂/(λ⁶*T³))*[(C₂/(λ*T))*(2*eᵃ¹/(eᵃ¹-1)-1)-2]
d²Ibb/dT² = C₁*(eᵃ¹/(eᵃ¹-1)²)*(a₁/(λ⁵*T²))*[a₁*(2*eᵃ¹/(eᵃ¹-1) -1)-2]
a₁=C₂/(λ*T)
a₂ = 1/(eᵃ¹-1)   #  1/expm1(a1)
a₃ = eᵃ¹/(eᵃ¹-1) #  exp(a)/expm1(a)
d²Ibb/dT² = C₁*a₂*a₃*(a₁/(λ⁵*T²))*[a₁*(2*a₃ - 1))-2]

Input :
            h  - to be filled, [Nx0]
            λ - wavelength in μm, [Nx0]
            T - tmperature in Kelvins
source
PlanckFunctions.∇²ₜibb!Method
∇²ₜibb!(h::AbstractVector{Float64},λ::AbstractVector{Float64},T::Float64,amat::AbstractMatrix{Float64})::Nothing

In-place bb intensity second order derivative with respect to temperature with 
intermediate matrix provided externally

d²Ibb/dT² = C₁*(eᵃ¹/(eᵃ¹-1)²)*(C₂/(λ⁶*T³))*[(C₂/(λ*T))*(2*eᵃ¹/(eᵃ¹-1)-1)-2]
d²Ibb/dT² = C₁*(eᵃ¹/(eᵃ¹-1)²)*(a₁/(λ⁵*T²))*[a₁*(2*eᵃ¹/(eᵃ¹-1) -1)-2]
a₁=C₂/(λ*T)
a₂ = 1/(eᵃ¹-1)   #  1/expm1(a1)
a₃ = eᵃ¹/(eᵃ¹-1) #  exp(a)/expm1(a)
d²Ibb/dT² = C₁*a₂*a₃*(a₁/(λ⁵*T²))*[a₁*(2*a₃ - 1))-2]

Input :
        h  - to be filled, [Nx0]
        λ - wavelength in μm, [Nx0]
        T - temperature in Kelvins
        amat - matrix of intermediate coefficients,  [Nx3]
source
PlanckFunctions.∇²ₜibb!Method
∇²ₜibb!(h::AbstractVector{Float64},λ::AbstractVector{Float64},T::Float64)

In-place bb intensity second order derivative with respect to temperature

d²Ibb/dT² = C₁*(eᵃ¹/(eᵃ¹-1)²)*(C₂/(λ⁶*T³))*[(C₂/(λ*T))*(2*eᵃ¹/(eᵃ¹-1)-1)-2]
d²Ibb/dT² = C₁*(eᵃ¹/(eᵃ¹-1)²)*(a₁/(λ⁵*T²))*[a₁*(2*eᵃ¹/(eᵃ¹-1) -1)-2]
a₁=C₂/(λ*T)
a₂ = 1/(eᵃ¹-1)   #  1/expm1(a1)
a₃ = eᵃ¹/(eᵃ¹-1) #  exp(a)/expm1(a)
d²Ibb/dT² = C₁*a₂*a₃*(a₁/(λ⁵*T²))*[a₁*(2*a₃ - 1))-2]

Input :
            h  - to be filled, [Nx0]
            λ - wavelength in μm, [Nx0]
            T - tmperature in Kelvins
source
PlanckFunctions.∇²ₜibb!Method
∇²ₜibb!(h::AbstractVector{Float64},T::Float64,amat::AbstractMatrix{Float64},∇i::AbstractVector{Float64})::Nothing

In-place bb intensity second order derivative with respect to temperature 
with provided both the intermediate matrix amat and the the Planck function first derivative

d²Ibb/dT² = C₁*(eᵃ¹/(eᵃ¹-1)²)*(C₂/(λ⁶*T³))*[(C₂/(λ*T))*(2*eᵃ¹/(eᵃ¹-1)-1)-2]
d²Ibb/dT² = C₁*(eᵃ¹/(eᵃ¹-1)²)*(a₁/(λ⁵*T²))*[a₁*(2*eᵃ¹/(eᵃ¹-1) -1)-2]
a₁=C₂/(λ*T)
a₂ = 1/(eᵃ¹-1)   #  1/expm1(a1)
a₃ = eᵃ¹/(eᵃ¹-1) #  exp(a)/expm1(a)
d²Ibb/dT² = C₁*a₂*a₃*(a₁/(λ⁵*T²))*[a₁*(2*a₃ - 1))-2]
    as far as 
        Ibb = (λ⁻⁵)* C₁*a₂
    and 
        dIbb/dT = C₁*a₃*a₂*a₁*(1/(λ⁵*T)) = a₃*a₁*Ibb/T 
    hense
        d²Ibb/dT² = C₁*a₂*a₃*a₁*(1/(λ⁵*T²))*[a₁*(2*a₃ - 1))-2] 
            = [a₃*a₁*Ibb/T^2]*[a₁*(2*a₃ - 1))-2] 
                = [(dIbb/dT)/T]*[a₁*(2*a₃ - 1))-2] 
Input :
        h  - to be filled, [Nx0]
        λ - wavelength in μm, [Nx0]
        amat - matrix of intermediate coefficients,  [Nx3]
        ∇i - vector of bb intensity first derivatives, [Nx0]
source
PlanckFunctions.∇²ₜibbMethod
∇²ₜibb(λ,T)

BB intensity second derivative with respect to temperature

d²Ibb/dT² = C₁*(eᵃ¹/(eᵃ¹-1)²)*(C₂/(λ⁶*T³))*[(C₂/(λ*T))*(2*eᵃ¹/(eᵃ¹-1)-1)-2]
d²Ibb/dT² = C₁*(eᵃ¹/(eᵃ¹-1)²)*(a₁/(λ⁵*T²))*[a₁*(2*eᵃ¹/(eᵃ¹-1) -1)-2]
a₁=C₂/(λ*T)
a₂ = 1/(eᵃ¹-1)   #  1/expm1(a1)
a₃ = eᵃ¹/(eᵃ¹-1) #  exp(a)/expm1(a)
d²Ibb/dT² = C₁*a₂*a₃*(a₁/(λ⁵*T²))*[a₁*(2*a₃ - 1))-2]

Input :
            λ - wavelength in μm
            T - tmperature in Kelvins
source
PlanckFunctions.∇ₗibbMethod
∇ₗibb(λ,T)


BB intensity first derivative with respect to the wavelength

Input:
    λ - wavelength, μm
    T - temperature, K
source
PlanckFunctions.∇ₜibb!Method
∇ₜibb!(g::AbstractMatrix,λ::AbstractVector,T::AbstractVector)

In-place BB intensity first derivative with respect to temperature
a₁=C₂/(λ*T)
a₂ = 1/(eᵃ¹-1)   #  1/expm1(a1)
a₃ = eᵃ¹/(eᵃ¹-1) #  exp(a)/expm1(a)
dIbb/dT = C₁*a₃*a₂*a₁*(1/(λ⁵*T))
Input:
    g - vector to be filled, [Nx0]
    λ - wavelength in μm, [Nx0]
    T - temperature in Kelvins
source
PlanckFunctions.∇ₜibb!Method
∇ₜibb!(g::AbstractVector,λ::AbstractVector,T,amat::AbstractMatrix)

In-place bb intensity first derivative with respect to temperature
with externally provided amat  - matrix with columns a₁,a₂,a₃

dIbb/dT = C₁*(eᵃ¹/(eᵃ¹-1)²)*(C₂/(λ⁶*T²))
a₁=C₂/(λ*T)
a₂ = 1/(eᵃ¹-1)   #  1/expm1(a1)
a₃ = eᵃ¹/(eᵃ¹-1) #  exp(a)/expm1(a)
dIbb/dT = C₁*a₃*a₂*a₁*(1/(λ⁵*T))

Input:
    g - to be filled, [Nx0]
    λ - wavelength in μm, [Nx0]
    T - temperature in Kelvins
    amat - matrix of intermediate coefficients, [Nx3]
source
PlanckFunctions.∇ₜibb!Method
∇ₜibb!(g::AbstractVector,T, amat::AbstractMatrix,i::AbstractVector)::Nothing

In-place bb intensity first derivative with respect to temperature
with externally provided amat  - matrix with columns a₁,a₂,a₃

dIbb/dT = C₁*(eᵃ¹/(eᵃ¹-1)²)*(C₂/(λ⁶*T²))
a₁=C₂/(λ*T)
a₂ = 1/(eᵃ¹-1)   #  1/expm1(a1)
a₃ = eᵃ¹/(eᵃ¹-1) #  exp(a)/expm1(a)
dIbb/dT = C₁*a₃*a₂*a₁*(1/(λ⁵*T))
as far as Ibb = C₁*a₂/λ⁵
dIbb/dT = a₃*a₁*C₁*(a₂/λ⁵)*(1/T)=a₃*a₁*Ibb/T

Input:
    g - to be filled, [Nx0]
    λ - wavelength in μm, [Nx0]
    T - temperature in Kelvins
    amat - matrix of intermediate coefficients, [Nx3]
source
PlanckFunctions.∇ₜibbMethod
∇ₜibb(λ::AbstractVector,T,amat::AbstractMatrix)

BB intensity first derivative with respect to temperature
with externally provided matrix of intermediate coefficients
dIbb/dT = C₁*(eᵃ¹/(eᵃ¹-1)²)*(C₂/(λ⁶*T²))
a₁=C₂/(λ*T)
a₂ = 1/(eᵃ¹-1)   
a₃ = eᵃ¹/(eᵃ¹-1) 
dIbb/dT = C₁*a₃*a₂*a₁*(1/(λ⁵*T))
Input:
    λ - wavelength in μm, [Nx0]
    T - temperature in Kelvins
    amat - matrix of intermediate coefficients, [Nx3]
source
PlanckFunctions.∇ₜibbMethod
∇ₜibb(λ,T)

BB intensity first derivative with respect to temperature
dIbb/dT = C₁*(eᵃ¹/(eᵃ¹-1)²)*(C₂/(λ⁶*T²))
a₁=C₂/(λ*T)
a₂ = 1/(eᵃ¹-1)   #  1/expm1(a1)
a₃ = eᵃ¹/(eᵃ¹-1) #  exp(a)/expm1(a)
dIbb/dT = C₁*a₃*a₂*a₁*(1/(λ⁵*T))
Input:
    λ - wavelength in μm
    T - temperature in Kelvins
source
PlanckFunctions.∫ibbₗMethod
∫ibbₗ(T;λₗ=0.0,λᵣ=Inf,tol=1e-6)

Relative (with respect to the integral power in the whole spectrum)
integral intensity of bb in the spectral range λₗ...λᵣ (by default the range is 0...inf)

Input:
    T - temperature,Kelvins
    (optional)
    λₗ - left wavelength boundary, μm
    λᵣ - right wavelength boundary, μm
    tol - intergation tolerance
source