PolarizedBRF.jl

PolarizedBRF.jl is a Julia package for calculating polarized bidrectional reflectance factors. It is largely based on M. Mishchenko et al.'s code.

PolarizedBRF.F_from_expansionMethod

Calculate scattering matrix for given angles using the expansion coefficients.

F_from_expansion(coeff, θ; output_dataframe, data)
  • coeff: Expansion coefficients, which is an lmax1 x 6 matrix.
  • θ: The angles in ascending order at which the scattering matrix is calculated.

Optional parameters:

  • output_dataframe: Whether or not to output the scattering matrix in the DataFrame format. Default is false, and an Nθ x 6 matrix will be outputed.
  • data: Use precalculated Wigner d data. By default new data will be calculated.
PolarizedBRF.evaluateMethod

Evaluate the refletion matrix $\mathbf{R}(\mu,\phi;\mu_0,\phi_0)$ using built interpolations.

evaluate(itps, μ, ϕ, μ₀, ϕ₀)
PolarizedBRF.expandMethod

Expand given scattering matrix (6-independent columns: a1, a2, a3, a4, b1, b2) to General Spherical Function (GSF) coefficients.

expand(F, θ₀; smax, smin, ngauss, double_threshold, adding_step, stop_threshold, cutoff, interpolation_strategy, ε, θₑ, error_type, _kwargs...)
  • F: The scattering matrix columns, which should be an Nθ x 6 matrix.
  • θ₀: The scattering angles. If the maximum value is larger than π, then the angles are assumed to be in degrees, otherwise the angles are assumed to be in radians.

Optional parameters:

  • smax: Highest level for expansion. Default is -1, and in such case, the highest level is determined via trial-and-error until the desired precision (defined by ε and error_type) is reached.
  • smin: The start point of iteration. Default is 30.
  • ngauss: The function for determine the number of Gaussian quadrature points. It takes the current expansion level plus 1, s + 1, as the input. Default is identity.
  • double_threshold: The highest level before which s is doubled each time. Default is 5000.
  • adding_step: The step size in the adding phase. Default is 100.
  • stop_threshold: The highest permitted expansion level. Default is 20000.
  • interpolation_strategy: The interpolation strategy. Default is PolarizedBRF.Linear, which uses LinearInterpolation from Interpolations.jl. Optional is PolarizedBRF.Spline, which uses Spline1D from Dierckx.jl.
  • cutoff: Cut off the expansion coefficients when the absolute value of α₁ at this level is smaller than cutoff. Default is 0.0, meaning no cutoff will be applied. Note that Ito et al. (2018) used cutoff = 1e-8.
  • ε: The desired error bound. Default is 0.02.
  • θₑ: The angles used for error checking. Default is θ₀.
  • error_type: The error type. Default is PolarizedBRF.AbsoluteError, other options include PolarizedBRF.RelativeError and PolarizedBRF.RootMeanSquaredError.
PolarizedBRF.get_itpMethod

Build interpolations from the calculated Fourier coefficients.

get_itp(x, R)
PolarizedBRF.ssfMethod

Calculate the static structure factor.

ssf(f, r, θ)
  • f: Volume fraction.
  • r: Volumetrically equivalent radius.
  • θ: Scattering angle in radians.
PolarizedBRF.ssf_correctionMethod

Apply the static structure factor correction to the given expansion coeffcients.

ssf_correction(coeff, Cext, Csca; volume_fraction, r, Nθ, Nm, equidistant, output_dataframe, strategy, kwargs...)
  • coeff: The original expansion coeffcients, an lmax1 x 6 matrix.
  • Cext: Extinction cross section.
  • Csca: Scattering cross section.

Optional parameters:

  • volume_fraction: The volume fraction taken by the scatterers. Default is 0.2.
  • r: The volumetrically equivalent/effective radius of the scatterer.
  • : Number of angles required for the output scattering matrix. Default is 181.
  • Nm: Number of angles used for the intermediate step. Default is 2lmax.
  • equidistant: Whether or not the output angles should be equidistant. Default is true. When set to false, Gaussian quadrature nodes (for cosθ) will be used instead.
  • output_dataframe: Whether or not to output the scattering matrix in the DataFrame format. Default is false, and an Nθ x 6 matrix will be outputed.
  • strategy: Special strategy for expansion. Default is PolarizedBRF.StandardExpansion, and the expansion can be controlled via the keyword arguments specified for expand(). If use PolarizedBRF.Ito18, then the strategy used in Ito et al. (2018) will be used instead, in which the highest expansion level is set to 4 * lmax and the expansion coefficients are cut off at the threshold 1e-8.
  • All keyword arguments for expand() also apply here.
PolarizedBRF.Wrapper.run_pbrfMethod

Calculate the Fourier coefficients of the refletion matrix using PolarizedBRF.

run_pbrf(ω, ngauss, coeff; ε, x, w, mode)
  • ω is the single scattering albedo, which should be within (0, 1].
  • ngauss determines the number of points used in integration.
  • coeff is the expansion coefficients, which should be an R x 6 matrix.

Optional parameters:

  • ε is the threshold for convergence check. Detault is 1e-7.
  • mode determines the way to do the quadrature. Default is PolarizedBRF.StandardQuadrature (NQUADR=2 in the original code), which is a normal Gaussian-Legendre quadrature within [0, 1]. Other options are PolarizedBRF.NormalOrientedQuadrature (NQUADR=1 in the original code), PolarizedBRF.ExplicitNormalQuadrature (NQUADR=3 in the original code) and PolarizedBRF.CustomQuadrature. Note that:
    • PolarizedBRF.StandardQuadrature might not be suitable if the refletion in the normal direction is needed.
    • For PolarizedBRF.CustomQuadrature, you need to input suitable nodes x and weights w yourself.
  • x is the custom quadrature nodes within [0, 1]. Only used when mode is PolarizedBRF.CustomQuadrature.
  • w is the corresponding quadrature weights. Only used when mode is PolarizedBRF.CustomQuadrature.

Results:

  • R is a 4 x 4 x NG x NG x LMAX1 array, storing all the Fourier coefficients.
  • x is the quadrature nodes.
  • w is the quadrature weights.