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_expansionPolarizedBRF.Wrapper.ddPolarizedBRF.Wrapper.run_pbrfPolarizedBRF.evaluatePolarizedBRF.expandPolarizedBRF.get_itpPolarizedBRF.ssfPolarizedBRF.ssf_correctionPolarizedBRF.ssf_hunterPolarizedBRF.ssf_sticky
PolarizedBRF.F_from_expansion — MethodCalculate scattering matrix for given angles using the expansion coefficients.
F_from_expansion(coeff, θ; output_dataframe, data)
coeff: Expansion coefficients, which is anlmax1 x 6matrix.θ: 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 theDataFrameformat. Default isfalse, and anNθ x 6matrix will be outputed.data: Use precalculated Wigner d data. By default new data will be calculated.
PolarizedBRF.evaluate — MethodEvaluate the refletion matrix $\mathbf{R}(\mu,\phi;\mu_0,\phi_0)$ using built interpolations.
evaluate(itps, μ, ϕ, μ₀, ϕ₀)
PolarizedBRF.expand — MethodExpand 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 anNθ x 6matrix.θ₀: 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εanderror_type) is reached.smin: The start point of iteration. Default is30.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 isidentity.double_threshold: The highest level before whichsis doubled each time. Default is5000.adding_step: The step size in the adding phase. Default is100.stop_threshold: The highest permitted expansion level. Default is20000.interpolation_strategy: The interpolation strategy. Default isPolarizedBRF.Linear, which usesLinearInterpolationfromInterpolations.jl. Optional isPolarizedBRF.Spline, which usesSpline1DfromDierckx.jl.cutoff: Cut off the expansion coefficients when the absolute value ofα₁at this level is smaller thancutoff. Default is0.0, meaning no cutoff will be applied. Note that Ito et al. (2018) usedcutoff = 1e-8.ε: The desired error bound. Default is0.02.θₑ: The angles used for error checking. Default isθ₀.error_type: The error type. Default isPolarizedBRF.AbsoluteError, other options includePolarizedBRF.RelativeErrorandPolarizedBRF.RootMeanSquaredError.
PolarizedBRF.get_itp — MethodBuild interpolations from the calculated Fourier coefficients.
get_itp(x, R)
PolarizedBRF.ssf — MethodCalculate the static structure factor.
ssf(f, r, θ)
f: Volume fraction.r: Volumetrically equivalent radius.θ: Scattering angle in radians.
PolarizedBRF.ssf_correction — MethodApply 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(lmax + 1) x 6matrix.Cext: Extinction cross section.Csca: Scattering cross section.
Optional parameters:
volume_fraction: The volume fraction taken by the scatterers. Default is0.2.r: The volumetrically equivalent/effective radius of the scatterer.τ: The stickiness factor. Default is0, meaning stickiness is not considered. The valid range isτ ≥ (2 - √2) / 6.Nθ: Number of angles required for the output scattering matrix. Default is181.Nm: Number of angles used for the intermediate step. Default is2lmax.equidistant: Whether or not the output angles should be equidistant. Default istrue. When set tofalse, Gaussian quadrature nodes (for cosθ) will be used instead.output_dataframe: Whether or not to output the scattering matrix in theDataFrameformat. Default isfalse, and anNθ x 6matrix will be outputed.strategy: Special strategy for expansion. Default isPolarizedBRF.StandardExpansion, and the expansion can be controlled via the keyword arguments specified forexpand(). If usePolarizedBRF.Ito18, then the strategy used in Ito et al. (2018) will be used instead, in which the highest expansion level is set to4 * lmaxand the expansion coefficients are cut off at the threshold1e-8.- All keyword arguments for
expand()also apply here.
PolarizedBRF.ssf_hunter — MethodCalculate the static structure factor using Hunter (2001)'s formula, which is almost equivalent to ssf(), but does not treat $|4r\sin\frac{\theta}{2}|<0.01$ cases specially.
ssf_hunter(f, r, θ)f: Volume fraction.r: Volumetrically equivalent radius.θ: Scattering angle in radians.
PolarizedBRF.ssf_sticky — MethodCalculate SSF for sticky particles.
ssf_sticky(f, r, θ, τ)
See Eqs. (8.4.20–22) in:
Tsang, L., Kong, J.A., Ding, K.-H., Ao, C.O., 2001. Scattering of Electromagnetic Waves: Numerical Simulations. John Wiley & Sons, Inc., New York, USA. https://doi.org/10.1002/0471224308
f: Volume fraction.r: Volumetrically equivalent radius.θ: Scattering angle in radians.τ: Stickiness parameter.
PolarizedBRF.Wrapper.dd — MethodCalculate $d_{m0}^s(x)$, $d_{m,-2}^s(x)$ and $d_{m,2}^s(x)$ for $s\in[0,l_{\max}]$.
dd(x, lmax, m)
PolarizedBRF.Wrapper.run_pbrf — MethodCalculate 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].ngaussdetermines the number of points used in integration.coeffis the expansion coefficients, which should be anR x 6matrix.
Optional parameters:
εis the threshold for convergence check. Detault is1e-7.modedetermines the way to do the quadrature. Default isPolarizedBRF.StandardQuadrature(NQUADR=2in the original code), which is a normal Gaussian-Legendre quadrature within[0, 1]. Other options arePolarizedBRF.NormalOrientedQuadrature(NQUADR=1in the original code),PolarizedBRF.ExplicitNormalQuadrature(NQUADR=3in the original code) andPolarizedBRF.CustomQuadrature. Note that:PolarizedBRF.StandardQuadraturemight not be suitable if the refletion in the normal direction is needed.- For
PolarizedBRF.CustomQuadrature, you need to input suitable nodesxand weightswyourself.
xis the custom quadrature nodes within[0, 1]. Only used whenmodeisPolarizedBRF.CustomQuadrature.wis the corresponding quadrature weights. Only used whenmodeisPolarizedBRF.CustomQuadrature.
Results:
Ris a4 x 4 x NG x NG x LMAX1array, storing all the Fourier coefficients.xis the quadrature nodes.wis the quadrature weights.