pytrunc.utils module#

pytrunc.utils.bessel_j1(x, acc=1e-08, max_iter=50)[source]#

The Bessel first kind function J1(x) of order 1

Parameters:
x1-D ndarray

The variable x of J1(x)

accfloat, optional

The tolerance for numerical errors. Default is 1e-8.

max_iterfloat, optional

The maximun number of iteration trying to improve the error accuracy

Returns:
J11-D ndarray

The values of the Bessel function J1(x)

Notes

The scipy equivalent -> scipy.special.j1(x)

pytrunc.utils.bessel_j1_derivative(x, acc=1e-08, max_iter=50)[source]#

Compute the Bessel first kind derivative of order 1 d(J1(x))

Parameters:
x1-D ndarray

The variable x of d(J1(x))

accfloat, optional

The tolerance for numerical errors. Default is 1e-8.

max_iterfloat, optional

The maximun number of iteration trying to improve the error accuracy

Returns:
dj11-D ndarray

The values of the Bessel function derivative d(J1(x))

Notes

The scipy equivalent -> scipy.special.jvp(1,x)

pytrunc.utils.bessel_j1_roots(nb_roots, acc=1e-08, max_iter=50)[source]#

Find roots of Bessel first kind function j1(x) using Newton-Raphson iteration

  • First k approximations equation ji_roots ~ π * ( k + π/2 + π/4 ). See ref[1]

Parameters:
nb_rootsint

The number of j1(x)=0 to find

accfloat, optional

The tolerance for numerical errors. Default is 1e-8.

max_iterfloat, optional

The maximun number of iteration trying to improve the error accuracy

Returns:
roots1-D ndarray

The roots of the function f(x)

Notes

The scipy equivalent -> scipy.special.jn_zeros(1,x)

References

  • [1] Baricz, Á., Kumar, P., & Ponnusamy, S. (2025). Asymptotic behavior of zeros

    of Bessel function derivatives. arXiv preprint arXiv:2510.12353.

pytrunc.utils.bessel_jn(x, n, acc=1e-08, max_iter=50)[source]#

The Bessel first kind function Jn(x) of order n

Parameters:
x1-D ndarray

The variable x of Jn(x)

nfloat

The Bessel first kind function order

accfloat, optional

The tolerance for numerical errors. Default is 1e-8.

max_iterfloat, optional

The maximun number of iteration trying to improve the error accuracy

Returns:
Jn1-D ndarray

The values of the Bessel function Jn(x)

Notes

The scipy equivalent -> scipy.special.jn(n,x)

pytrunc.utils.integrate_lobatto(f, x, lp=None, xk=None, wk=None, assume_sorted=False)[source]#

Integrate using lobatto quadrature

Parameters:
f1-D ndarray

The ordinates of the function (array to be integrated).

x1-D ndarray

The abscissas

lpNone | int, optional

The number of lobatto points for the integration. If None lp = len(x).

xkNone | 1-D ndarray

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

wkNone | 1-D ndarray

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

assume_sortedbool, optional

If True, the x array is assumed to be sorted in ascending order.

Return
——
intfloat

The estimated integral calculated using the Lobatto quadrature

pytrunc.utils.legendre_polynomials(n, x)[source]#

Use the recursion formulas to compute the Legendre polynomials Pn(x)

  • see Eq. 9 in [1]

Parameters:
nint

The Legendre polynomial order

x1-D ndarray

The x values of Pn(x)

Returns:
P1-D ndarray

The Legendre series

Notes

The numpy equivalent -> numpy.polynomial.legendre.Legendre.basis(n)(x)

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.utils.legendre_polynomials_derivative(n, x)[source]#

Use the recursion formulas to compute the derivative Legendre polynomials d(Pn(x))

  • see Eq. 10 in [1]

Parameters:
nint

The Legendre polynomial order

x1-D ndarray

The x values of d(Pn(x))

Returns:
P1-D ndarray

The derivative Legendre series

Notes

The numpy equivalent -> numpy.polynomial.legendre.Legendre.basis(n).deriv(1)(x)

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.utils.legendre_polynomials_derivative_roots(n, acc=1e-08, max_iter=50)[source]#

Find roots of legendre polynomial derivative d(Pn(x)), for x > -1 and x < 1

  • Use of Newton-raphson iteration as [1], see Eq.7 and 8.

Parameters:
nint

The Legendre polynomial order

accfloat, optional

The tolerance for numerical errors. Default is 1e-8.

max_iterfloat, optional

The maximun number of iteration trying to improve the error accuracy

Returns:
roots1-D ndarray

The roots of legendre polynomial derivative d(Pn(x))

Notes

  • The numpy equivalent -> Legendre.basis(n).deriv().roots()

  • Faster than the numpy equivalent!

References

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

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

pytrunc.utils.legendre_polynomials_second_derivative(n, x)[source]#

Use the recursion formulas to compute the second derivative Legendre polynomials d²(Pn(x))

  • see Eq. 11 in [1]

Parameters:
nint

The Legendre polynomial order

x1-D ndarray

The x values of d²(Pn(x))

Returns:
P1-D ndarray

The second derivative Legendre series

Notes

The numpy equivalent -> numpy.polynomial.legendre.Legendre.basis(n).deriv(2)(x)

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.utils.quadrature_lobatto(abscissa_min=-1, abscissa_max=1, n=100)[source]#

Compute the abscissas (sample points) and weigths for Lobatto quadrature.

Parameters:
abscissa_minfloat

The min fixed abscissa

abscissa_maxfloat

The max fixed abscissa

nint, optional

The number of abscissas / weights

Returns:
abscissas1-D ndarray

The Lobatto quadrature abscissas

weights1-D ndarray

The Lobatto quadrature weights

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.