pyrt.henyey_greenstein_legendre_coefficients#
- pyrt.henyey_greenstein_legendre_coefficients(asymmetry_parameter, n_coefficients)[source]#
Get the Legendre coefficients of a Henyey-Greenstein phase function.
- Parameters:
asymmetry_parameter (float) – The Henyey-Greenstein asymmetry parameter. Must be between -1 and 1.
n_coefficients (int) – The number of coefficients to keep.
- Returns:
1-dimensional arrray of Legendre coefficients up to the specificed number of coefficients.
- Return type:
np.ndarray
Notes
The Henyey-Greenstein phase function can be decomposed as follows:
\[p(\mu) = \sum_{n=0}^{\infty} (2n + 1)g^n P_n(\mu)\]where \(p\) is the phase function, \(\mu\) is the cosine of the scattering angle, \(n\) is the moment number, \(g\) is the asymmetry parameter, and \(P_n(\mu)\) is the \(n\)th Legendre polynomial.
Examples
Get the first 129 coefficients of the Henyey-Greenstein phase function for a given asymmetry parameter.
>>> import numpy as np >>> import pyrt >>> g = 0.5 >>> coeff = pyrt.henyey_greenstein_legendre_coefficients(g, 129) >>> coeff.shape (129,)
Construct a Henyey-Greenstein phase function, decompose it, and see how this result compares to the analytic decomposition performed above.
>>> ang = np.linspace(0, 180, num=181) >>> pf = pyrt.construct_henyey_greenstein(g, ang) * 4 * np.pi # normalize it >>> lc = pyrt.decompose(pf, ang, 129) >>> print(np.amax(np.abs(lc - coeff)) < 1e-10) True