Functions for blackbody thermal emission (and its derivatives) calculation

Main.BandPyrometry.Planck.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
Main.BandPyrometry.Planck.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
Main.BandPyrometry.Planck.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
Main.BandPyrometry.Planck.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
Main.BandPyrometry.Planck.band_powerMethod
band_power(T;λₗ=0.0,λᵣ=Inf,tol=1e-6)

Total bb with temperature T integral intensity within 
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
Main.BandPyrometry.Planck.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
Main.BandPyrometry.Planck.ibb!Method
ibb!(i::AbstractVector,λ::AbstractVector,T::Float64)

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

Blackbody intensity with intermediate matrix provided externally, [W/m2-sr-mkm]
    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
Main.BandPyrometry.Planck.ibbMethod
ibb(λ::AbstractVector,T::AbstractVector)

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

Blackbody intensity , [W/m²⋅sr⋅μm]
Ibb = (λ⁻⁵)* C₁/(eᵃ¹-1) , where a₁=C₂/(λ*T)
Input:
    λ - wavelength in μm
    T - temperature in Kelvins
source
Main.BandPyrometry.Planck.∇²ₜ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
Main.BandPyrometry.Planck.∇²ₜ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
Main.BandPyrometry.Planck.∇²ₜ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
Main.BandPyrometry.Planck.∇²ₜ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
Main.BandPyrometry.Planck.∇²ₜ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
Main.BandPyrometry.Planck.∇ₜ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
Main.BandPyrometry.Planck.∇ₜ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
Main.BandPyrometry.Planck.∇ₜ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
Main.BandPyrometry.Planck.∇ₜ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
Main.BandPyrometry.Planck.∇ₜ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
Main.BandPyrometry.Planck.∫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