MieScattering.MieScattering
MieScattering.D_calc
MieScattering.D_downwards!
MieScattering.D_upwards!
MieScattering.Lentz_Dn
MieScattering.ez_intensities
MieScattering.ez_mie
MieScattering.generate_mie_costheta
MieScattering.i_par
MieScattering.i_per
MieScattering.i_unpolarized
MieScattering.mie
MieScattering.mie_An_Bn
MieScattering.mie_S1_S2
MieScattering.mie_cdf
MieScattering.mie_mu_with_uniform_cdf
MieScattering.mie_phase_matrix
MieScattering.normalization_factor
MieScattering.small_conducting_mie
MieScattering.small_mie
MieScattering.small_mie_S1_S2
MieScattering.small_mie_conducting_S1_S2
MieScattering.MieScattering
— ModuleMie scattering calculations for perfect spheres based on miepython. Extensive documentation is at (https://miepython.readthedocs.io).
MieScattering.jl
is a Julia package to calculate light scattering of a plane wave by non-np.absorbing, partially-np.absorbing, or perfectly conducting spheres.
The extinction efficiency, scattering efficiency, backscattering, and scattering asymmetry for a sphere with complex index of refraction m, diameter d, and wavelength lambda can be found by:
qext, qsca, qback, g = ez_mie(m, d, λ0)
The normalized scattering values for angles µ=cos(θ) are:
Ipar, Iper = ez_intensities(m, d, λ0, µ)
If the size parameter is known, then use:
mie(m, x)
Mie scattering amplitudes S1 and S2 (complex numbers):
mie_S1_S2(m, x, μ)
Normalized Mie scattering intensities for angles µ=cos(θ):
i_per(m, x, µ)
i_par(m, x, µ)
i_unpolarized(m, x, µ)
MieScattering.D_calc
— MethodD_calc(m, x, N)
Compute the logarithmic derivative using best method.
Parameters
m
: the complex index of refraction of the spherex
: the size parameter of the sphereN
: order of Ricatti-Bessel function
Output
The values of the Ricatti-Bessel function for orders from 0 to N.
MieScattering.D_downwards!
— MethodD_downwards!(z, N, D)
Compute the logarithmic derivative by downwards recurrence.
Parameters
z
: function argumentN
: order of Ricatti-Bessel functionD
: gets filled with the Ricatti-Bessel function values for orders from 0 to N for an argument z using the downwards recurrence relations.
MieScattering.D_upwards!
— MethodD_upwards!(z, N, D)
Compute the logarithmic derivative by upwards recurrence.
Parameters
z
: function argumentN
: order of Ricatti-Bessel functionD
: gets filled with the Ricatti-Bessel function values for orders from 0 to N for an argument z using the upwards recurrence relations.
MieScattering.Lentz_Dn
— MethodLentz_Dn(z, N)
Compute the logarithmic derivative of the Ricatti-Bessel function.
Parameters
z
: function argumentN
: order of Ricatti-Bessel function
Output
This returns the Ricatti-Bessel function of order N with argument z using the continued fraction technique of Lentz, Appl. Opt., 15, 668-671, (1976).
MieScattering.ez_intensities
— Functionez_intensities(m, d, λ0, μ, n_env = 1.0, norm = :albedo)
Return the scattered intensities from a sphere. These are the scattered intensities in a plane that is parallel (ipar) and perpendicular (iper) to the field of the incident plane wave. The scattered intensity is normalized such that the integral of the unpolarized intensity over 4𝜋 steradians is equal to the single scattering albedo. The scattered intensity has units of inverse steradians [1/sr]. The unpolarized scattering is the average of the two scattered intensities.
Parameters
m
: the complex index of refraction of the sphere [-]d
: the diameter of the sphere [same units as lambda0]λ0
: wavelength in a vacuum [same units as d]µ
: the cos(θ) of each direction desired [-]n_env
: real index of medium around sphere, optional.
Output
ipar
, iper
: scattered intensity in parallel and perpendicular planes [1/sr]
MieScattering.ez_mie
— Functionez_mie(m, d, λ0, n_env = 1.0)
Calculate the efficiencies of a sphere.
Parameters
m
: the complex index of refraction of the sphere [-]d
: the diameter of the sphere [same units as lambda0]λ0
: wavelength in a vacuum [same units as d]n_env
: real index of medium around sphere, optional.use_threads
(optional): Flag whether to use threads (default: true)
Output
qext
: the total extinction efficiency [-]qsca
: the scattering efficiency [-]qback
: the backscatter efficiency [-]g
: the average cosine of the scattering phase function [-]
MieScattering.generate_mie_costheta
— Methodgenerate_mie_costheta(μ_cdf)
Generate a new scattering angle using a cdf. A uniformly spaced cumulative distribution function (CDF) is needed. New random angles are generated by selecting a random interval μ[i] to μ[i+1] and choosing an angle uniformly distributed over the interval.
Parameters:
μ_cdf
: a cumulative distribution function
Output:
- an array of random scattering angle cosines based on the CDF supplied.
MieScattering.i_par
— Methodi_par(m, x, μ; norm = :albedo)
Return the scattered intensity in a plane parallel to the incident light. This is the scattered intensity in a plane that is perpendicular to the field of the incident plane wave. The intensity is normalized such that the integral of the unpolarized intensity over 4π steradians is equal to the single scattering albedo.
Parameters
m
: the complex index of refraction of the spherex
: the size parameter of the sphereµ
: the angles, cos(theta), to calculate intensities
Output
The intensity at each angle in the array µ. Units [1/sr]
MieScattering.i_per
— Methodi_per(m, x, μ; norm = :albedo)
Return the scattered intensity in a plane normal to the incident light. This is the scattered intensity in a plane that is perpendicular to the field of the incident plane wave. The intensity is normalized such that the integral of the unpolarized intensity over 4π steradians is equal to the single scattering albedo.
Parameters
m
: the complex index of refraction of the spherex
: the size parameter of the sphereµ
: the angles, cos(theta), to calculate intensities
Output
The intensity at each angle in the array µ. Units [1/sr]
MieScattering.i_unpolarized
— Methodi_unpolarized(m, x, μ; norm = :albedo)
Return the unpolarized scattered intensity at specified angles. This is the average value for randomly polarized incident light. The intensity is normalized such that the integral of the unpolarized intensity over 4π steradians is equal to the single scattering albedo.
Parameters
m
: the complex index of refraction of the spherex
: the size parameter of the sphereµ
: the angles, cos(theta), to calculate intensities
Output
The intensity at each angle in the array µ. Units [1/sr]
MieScattering.mie
— Functionmie(m, x)
Calculate the efficiencies for a sphere where m or x may be vectors.
Parameters
m
: the complex index of refraction of the spherex
: the size parameter of the sphere
Output
qext
: the total extinction efficiencyqsca
: the scattering efficiencyqback
: the backscatter efficiencyg
: the average cosine of the scattering phase function
MieScattering.mie_An_Bn
— Methodmie_An_Bn(m, x)
Compute arrays of Mie coefficients A and B for a sphere. This estimates the size of the arrays based on Wiscombe's formula. The length of the arrays is chosen so that the error when the series are summed is around 1e-6.
Parameters
m
: the complex index of refraction of the spherex
: the size parameter of the sphere
Output
An
,Bn
: arrays of Mie coefficents
MieScattering.mie_S1_S2
— MethodCalculate the scattering amplitude functions for spheres. The amplitude functions have been normalized so that when integrated over all 4*π
solid angles, the integral will be qext*pi*x^2
. The units are weird, $sr^{-0.5}$.
Parameters
m
: the complex index of refraction of the spherex
: the size parameter of the sphereµ
: the angles, cos($θ$), to calculate scattering amplitudesnorm
(optional): The normalization. Must be one of :albedo (default), :one, :four_pi, :qext,
:qsca, :bohren or :wiscombe"
use_threads
(optional): Flag whether to use threads (default: true)
Output
S1
, S2
: the scattering amplitudes at each angle µ [$sr^{-0.5}$]
MieScattering.mie_cdf
— Methodmie_cdf(m, x, num; norm = :albedo)
Create a CDF for unpolarized scattering uniformly spaced in cos(θ). The CDF covers scattered (exit) angles ranging from 180 to 0 degrees. (The cosines are uniformly distributed over -1 to 1.) Because the angles are uniformly distributed in cos(theta), the scattering function is not sampled uniformly and therefore huge array sizes are needed to adequately sample highly anisotropic phase functions. Since this is a cumulative distribution function, the maximum value should be 1.
Parameters
m
: the complex index of refraction of the spherex
: the size parameter of the spherenum
: length of desired CDF array
Output
µ
: array of cosines of anglescdf
: array of cumulative distribution function values
MieScattering.mie_mu_with_uniform_cdf
— Methodmie_cdf(m, x, num; norm = :albedo)
Create a CDF for unpolarized scattering for uniform CDF. The CDF covers scattered (exit) angles ranging from 180 to 0 degrees. (The cosines are uniformly distributed over -1 to 1.) These angles mu correspond to uniform spacing of the cumulative distribution function for unpolarized Mie scattering where cdf[i] = i/(num-1). This is a brute force implementation that solves the problem by calculating the CDF at many points and then scanning to find the specific angles that correspond to uniform interval of the CDF. Since this is a cumulative distribution function, the maximum value should be 1.
Parameters
m
: the complex index of refraction of the spherex
: the size parameter of the spherenum
: length of desired CDF array
Output
µ
: array of cosines of angles (irregularly spaced)cdf
: array of cumulative distribution function values
MieScattering.mie_phase_matrix
— Functionmie_phase_matrix(m, x, μ; norm=:albedo)
Calculate the phase scattering matrix.
The units are $sr^{-1}$. The phase scattering matrix is computed from the scattering amplitude functions, according to equations 5.2.105-6 in K. N. Liou (2002) - An Introduction to Atmospheric Radiation, Second Edition.
Parameters
m
: the complex index of refraction of the spherex
: the size parameter of the sphereμ
: the angles, cos(theta), at which to calculate the phase scattering matrix
Output
p
: The phase scattering matrix [$sr^{-1}$]
MieScattering.normalization_factor
— Methodnormalization_factor(a, b, x; norm)
Figure out scattering function normalization.
Parameters
m
: complex index of refraction of spherex
: dimensionless sphere sizenorm
: symbol describing type of normalization
Output
scaling factor needed for scattering function
MieScattering.small_conducting_mie
— Methodsmall_conducting_mie(m, x)
Calculate the efficiencies for a small conducting spheres. Typically used for small conducting spheres where x < 0.1
and real(m) == 0
.
Parameters
x
: the size parameter of the sphere
Output
qext
: the total extinction efficiencyqsca
: the scattering efficiencyqback
: the backscatter efficiencyg
: the average cosine of the scattering phase function
MieScattering.small_mie
— Methodsmall_mie(m, x)
Calculate the efficiencies for a small sphere. Typically used for small spheres where x<0.1
Parameters
m
: the complex index of refraction of the spherex
: the size parameter of the sphere
Output
qext
: the total extinction efficiencyqsca
: the scattering efficiencyqback
: the backscatter efficiencyg
: the average cosine of the scattering phase function
MieScattering.small_mie_S1_S2
— MethodCalculate the scattering amplitude functions for small spheres (x<0.1
). The amplitude functions have been normalized so that when integrated over all 4*π
solid angles, the integral will be qext*pi*x^2
. The units are weird, $sr^{-0.5}$
Parameters
m
: the complex index of refraction of the spherex
: the size parameter of the sphereµ
: the angles, cos($θ$), to calculate scattering amplitudes
Output
S1
, S2
: the scattering amplitudes at each angle µ [$sr^{-0.5}$]
MieScattering.small_mie_conducting_S1_S2
— Methodsmall_mie_conducting_S1_S2(m, x, μ)
Calculate the scattering amplitudes for small conducting spheres. The spheres are small perfectly conducting (reflecting) spheres (x<0.1
). The amplitude functions have been normalized so that when integrated over all 4𝜋
solid angles, the integral will be qext(𝜋x²
). The units are weird, $sr^{-0.5}$.
Parameters
m
: the complex index of refraction of the spherex
: the size parameter of the sphereµ
: the angles, cos($θ$), to calculate scattering amplitudes
Output
S1
, S2
: the scattering amplitudes at each angle µ [$sr^{-0.5}$]