pytrunc.phase module#

pytrunc.phase.calc_hg_moments(g, m_max)[source]#

Compute exact Henyey-Greenstein phase moments

  • see Eq.8 in [1]

Parameters:
gfloat

The Henyey-Greenstein parameter g (measures the asymmetry of the phase matrix)

m_maxint

The maximum moment number to compute, i.e., compute m[0], …, m[m_max]

Returns:
m1-D ndarray

The phase moment of size m_max + 1

References

  • [1] Kattawar, G. W. (1975). A three-parameter analytic phase function for multiple

    scattering calculations. Journal of Quantitative Spectroscopy and Radiative Transfer, 15(9), 839-849.

pytrunc.phase.calc_moments(phase, theta, m_max, method='lobatto', theta_unit='deg', normalize=False, xk=None, wk=None, pl_costh=None)[source]#

Calculate the phase matrix moments until m_max moment

Parameters:
phase1-D ndarray

The phase matrix

theta1-D ndarray

The phase matrix angles. See parameter theta_unit

m_maxint

The maximum moment number to compute, i.e., compute m[0], …, m[m_max]

methodstr, optional

The method used to calculate the moments, choices: - ‘lobatto’ -> Default value, very effcient if “gauss kind” theta distribution - ‘simpson’ -> efficient if regular theta distribution (use scipy.integrate.simpson)

theta_unitstr, optional

The unit for theta angles: - ‘deg’ (default value) - ‘rad’

normalizebool, optional

If normalize = True -> normalize such that first moment exactly = 1

xkNone | 1-D ndarray, optional

Force the Lobatto quadrature abscissas. Considered if wk is also provided.

wkNone | 1-D ndarray, optional

Force the Lobatto weights. Considered if xk is also provided.

pl_costhNone | 2-D ndarray, optional

Force the Legendre polynomials values of cos(theta). The shape must be (m_max+1, len(theta))

Returns
——-
m1-D ndarray

The computed phase moment of size m_max + 1

Notes

  • See Eq.A2 in ref[2] for moment computation using Lobatto quadrature in [0,π]

  • For Lobatto quadrature abcissas and weights calculation see ref[2]

References

  • [1] Michels, H. (1963). Abscissas and weight coefficients for Lobatto quadrature.

    Mathematics of Computation, 17(83), 237-244.

  • [2] Wiscombe, W. J. (1977). The delta-M method: Rapid yet accurate radiative flux calculations

    for strongly asymmetric phase functions. Journal of Atmospheric Sciences, 34(9), 1408-1422.

pytrunc.phase.calc_tthg_moments(g1, g2, f, m_max)[source]#

Compute exact Two-term Henyey-Greenstein phase moments

  • see Eq.11 in [1]

Parameters:
g1float

The first H-G term parameter g (forward part)

g2float

The second H-G term parameter g (backward part)

ffloat

The fraction parameter between the two H-G terms

m_maxint

The maximum moment number to compute, i.e., compute m[0], …, m[m_max]

Returns:
m1-D ndarray

The phase moment of size m_max + 1

References

  • [1] Kattawar, G. W. (1975). A three-parameter analytic phase function for multiple

    scattering calculations. Journal of Quantitative Spectroscopy and Radiative Transfer, 15(9), 839-849.

pytrunc.phase.henyey_greenstein(theta, g, theta_unit='deg', normalize=None)[source]#

Compute the Henyey-Greenstein phase matrix

Parameters:
theta1-D ndarray

The phase matrix angles. See parameter theta_unit

gfloat

The Henyey-Greenstein parameter g (measures the asymmetry of the phase matrix)

theta_unitstr, optional

The unit for theta angles: - ‘deg’ (default value) - ‘rad’

normalizeNone | float, optional

The normalization value of the integral ∫F_HG(θ)dcosθ, where F_HG(θ) is the phase matrix. The scipy simpson function is used for the normalization.

Return
——
F_HG1-D ndarray

The phase matrix.

Notes

The Henyey-Greenstein equation:

  • \(F_HG(θ) = (1/π) * [ (1 - g**2) / (1 + g**2 - (2*g*cos(θ)))**(3/2) ]\)

  • By default the integral ∫F(θ)dcosθ = 1/2π. The integral value can be different due to a very low dicretization of θ and/or a high g value. The use of the parameter normalize can be useful to renormalize the phase function.

References

Henyey, L. G., & Greenstein, J. L. (1941). Diffuse radiation in the galaxy. Astrophysical Journal, vol. 93, p. 70-83 (1941)., 93, 70-83.

other reference

pytrunc.phase.two_term_henyey_greenstein(theta, g1, g2, f, theta_unit='deg', normalize=None)[source]#

Compute the two-term Henyey-Greenstein phase matrix

Parameters:
theta1-D ndarray

The phase matrix angles. See parameter theta_unit

g1float

The first H-G term parameter g (forward part)

g2float

The second H-G term parameter g (backward part)

ffloat

The fraction parameter between the two H-G terms (see notes)

theta_unitstr, optional

The unit for theta angles: - ‘deg’ (default value) - ‘rad’

normalizeNone | float, optional

The normalization value of the integral ∫F(θ)dcosθ, where F(θ) is the phase matrix. The scipy simpson function is used for the normalization.

Return
——
F_TTHG1-D ndarray

The phase matrix.

Notes

The two term Henyey-Greenstein equation:

  • \(F_TTHG(θ) = f*F_HG1(θ) + (1-f)*F_HG2(θ)\)

  • By default the integral ∫F_TTHG(θ)dcosθ = 1/2π. The integral value can be different due to a very low dicretization of θ and/or a high g value. The use of the parameter normalize can be useful to renormalize the phase function.

References

Irvine, W. M. (1965). Multiple scattering by large particles (No. NASA-CR-64638).