prysm.polynomials#

Common Routines#

Various polynomials of optics.

prysm.polynomials.sum_of_2d_modes(modes, weights)#

Compute a sum of 2D modes.

Parameters:
  • modes (iterable) – sequence of ndarray of shape (k, m, n); a list of length k with elements of shape (m,n) works

  • weights (numpy.ndarray) – weight of each mode

Returns:

ndarray of shape (m, n) that is the sum of modes as given

Return type:

numpy.ndarry

prysm.polynomials.sum_of_2d_modes_backprop(modes, databar)#

Gradient backpropagation through sum_of_2d_modes.

Parameters:
  • modes (iterable) – sequence of ndarray of shape (k, m, n); a list of length k with elements of shape (m,n) works

  • databar (numpy.ndarray) – partial gradient backpropated up to the return of sum_of_2d_modes

Returns:

cumulative gradient through to the weights vector given to sum_of_2d_modes

Return type:

numpy.ndarry

prysm.polynomials.hopkins(a, b, c, r, t, H)#

Hopkins’ aberration expansion.

This function uses the “W020” or “W131” like notation, with Wabc separating into the a, b, c arguments. To produce a sine term instead of cosine, make a the negative of the order. In other words, for W222S you would use hopkins(2, 2, 2, …) and for W222T you would use hopkins(-2, 2, 2, …).

Parameters:
  • a (int) – azimuthal order

  • b (int) – radial order

  • c (int) – order in field (“H-order”)

  • r (numpy.ndarray) – radial pupil coordinate

  • t (numpy.ndarray) – azimuthal pupil coordinate

  • H (numpy.ndarray) – field coordinate

Returns:

polynomial evaluated at this point

Return type:

numpy.ndarray

prysm.polynomials.lstsq(modes, data)#

Least-Squares fit of modes to data.

Parameters:
  • modes (iterable) – modes to fit; sequence of ndarray of shape (m, n)

  • data (numpy.ndarray) – data to fit, of shape (m, n) place NaN values in data for points to ignore

Returns:

fit coefficients

Return type:

numpy.ndarray

Zernikes#

Zernike polynomials.

prysm.polynomials.zernike.zernike_norm(n, m)#

Norm of a Zernike polynomial with n, m indexing.

prysm.polynomials.zernike.zero_separation(n)#

Zero separation in normalized r based on radial order n.

prysm.polynomials.zernike.zernike_nm(n, m, r, t, norm=True)#

Zernike polynomial of radial order n, azimuthal order m at point r, t.

Parameters:
  • n (int) – radial order

  • m (int) – azimuthal order

  • r (numpy.ndarray) – radial coordinates

  • t (numpy.ndarray) – azimuthal coordinates

  • norm (bool, optional) – if True, orthonormalize the result (unit RMS) else leave orthogonal (zero-to-peak = 1)

Returns:

zernike mode of order n,m at points r,t

Return type:

numpy.ndarray

prysm.polynomials.zernike.zernike_nm_sequence(nms, r, t, norm=True)#

Zernike polynomial of radial order n, azimuthal order m at point r, t.

Parameters:
  • nms (iterable of tuple of int,) – sequence of (n, m); looks like [(1,1), (3,1), …]

  • r (numpy.ndarray) – radial coordinates

  • t (numpy.ndarray) – azimuthal coordinates

  • norm (bool, optional) – if True, orthonormalize the result (unit RMS) else leave orthogonal (zero-to-peak = 1)

Returns:

shape (k, n, m), with k = len(nms)

Return type:

numpy.ndarray

prysm.polynomials.zernike.zernike_nm_der(n, m, r, t, norm=True)#

Derivatives of Zernike polynomial of radial order n, azimuthal order m, w.r.t r and t.

Parameters:
  • n (int) – radial order

  • m (int) – azimuthal order

  • r (numpy.ndarray) – radial coordinates

  • t (numpy.ndarray) – azimuthal coordinates

  • norm (bool, optional) – if True, orthonormalize the result (unit RMS) else leave orthogonal (zero-to-peak = 1)

Returns:

dZ/dr, dZ/dt

Return type:

numpy.ndarray, numpy.ndarray

prysm.polynomials.zernike.zernike_nm_der_sequence(nms, r, t, norm=True)#

Derivatives of Zernike polynomial of radial order n, azimuthal order m, w.r.t r and t.

Parameters:
  • nms (iterable) – sequence of [(n, m)] radial and azimuthal orders

  • m (int) – azimuthal order

  • r (numpy.ndarray) – radial coordinates

  • t (numpy.ndarray) – azimuthal coordinates

  • norm (bool, optional) – if True, orthonormalize the result (unit RMS) else leave orthogonal (zero-to-peak = 1)

Returns:

length (len(nms)) list of (dZ/dr, dZ/dt)

Return type:

list

prysm.polynomials.zernike.nm_to_fringe(n, m)#

Convert (n,m) two term index to Fringe index.

prysm.polynomials.zernike.nm_to_ansi_j(n, m)#

Convert (n,m) two term index to ANSI single term index.

prysm.polynomials.zernike.ansi_j_to_nm(idx)#

Convert ANSI single term to (n,m) two-term index.

prysm.polynomials.zernike.noll_to_nm(idx)#

Convert Noll Z to (n, m) two-term index.

prysm.polynomials.zernike.fringe_to_nm(idx)#

Convert Fringe Z to (n, m) two-term index.

prysm.polynomials.zernike.zernikes_to_magnitude_angle_nmkey(coefs)#

Convert Zernike polynomial set to a magnitude and phase representation.

Parameters:

coefs (list of tuples) – a list looking like[(1,2,3),] where (1,2) are the n, m indices and 3 the coefficient

Returns:

dict keyed by tuples of (n, |m|) with values of (rho, phi) where rho is the magnitudes, and phi the phase

Return type:

dict

prysm.polynomials.zernike.zernikes_to_magnitude_angle(coefs)#

Convert Zernike polynomial set to a magnitude and phase representation.

This function is identical to zernikes_to_magnitude_angle_nmkey, except its keys are strings instead of (n, |m|)

Parameters:

coefs (list of tuples) – a list looking like[(1,2,3),] where (1,2) are the n, m indices and 3 the coefficient

Returns:

dict keyed by friendly name strings with values of (rho, phi) where rho is the magnitudes, and phi the phase

Return type:

dict

prysm.polynomials.zernike.nm_to_name(n, m)#

Convert an (n,m) index into a human readable name.

Parameters:
  • n (int) – radial polynomial order

  • m (int) – azimuthal polynomial order

Returns:

a name, np.g. Piston or Primary Spherical

Return type:

str

prysm.polynomials.zernike.top_n(coefs, n=5)#

Identify the top n terms in the wavefront expansion.

Parameters:
  • coefs (dict) – keys of (n,m), values of magnitudes, e.g. {(3,1): 2} represents 2 of primary coma

  • n (int, optional) – identify the top n terms.

Returns:

list of tuples (magnitude, index, term)

Return type:

list

prysm.polynomials.zernike.barplot(coefs, names=None, orientation='h', buffer=1, zorder=3, number=True, offset=0, width=0.8, fig=None, ax=None)#

Create a barplot of coefficients and their names.

Parameters:
  • coefs (dict) – with keys of Zn, values of numbers

  • names (dict) – with keys of Zn, values of names (e.g. Primary Coma X)

  • orientation (str, {'h', 'v', 'horizontal', 'vertical'}) – orientation of the plot

  • buffer (float, optional) – buffer to use around the left and right (or top and bottom) bars

  • zorder (int, optional) – zorder of the bars. Use zorder > 3 to put bars in front of gridlines

  • number (bool, optional) – if True, plot numbers along the y=0 line showing indices

  • offset (float, optional) – offset to apply to bars, useful for before/after Zernike breakdowns

  • width (float, optional) – width of bars, useful for before/after Zernike breakdowns

  • fig (matplotlib.figurnp.Figure) – Figure containing the plot

  • ax (matplotlib.axes.Axis) – Axis containing the plot

Returns:

  • fig (matplotlib.figurnp.Figure) – Figure containing the plot

  • ax (matplotlib.axes.Axis) – Axis containing the plot

prysm.polynomials.zernike.barplot_magnitudes(magnitudes, orientation='h', sort=False, buffer=1, zorder=3, offset=0, width=0.8, fig=None, ax=None)#

Create a barplot of magnitudes of coefficient pairs and their names.

e.g., astigmatism will get one bar.

Parameters:
  • magnitudes (dict) – keys of names, values of magnitudes. E.g., {‘Primary Coma’: 1234567}

  • orientation (str, {'h', 'v', 'horizontal', 'vertical'}) – orientation of the plot

  • sort (bool, optional) – whether to sort the zernikes in descending order

  • buffer (float, optional) – buffer to use around the left and right (or top and bottom) bars

  • zorder (int, optional) – zorder of the bars. Use zorder > 3 to put bars in front of gridlines

  • offset (float, optional) – offset to apply to bars, useful for before/after Zernike breakdowns

  • width (float, optional) – width of bars, useful for before/after Zernike breakdowns

  • fig (matplotlib.figure.Figure) – Figure containing the plot

  • ax (matplotlib.axes.Axis) – Axis containing the plot

Returns:

  • fig (matplotlib.figure.Figure) – Figure containing the plot

  • ax (matplotlib.axes.Axis) – Axis containing the plot

XY#

XY polynomials.

prysm.polynomials.xy.j_to_xy(j)#

Convert a mono-index j into the x and y powers.

counts from j=2 for x^1 y^0 and j=3 for x^0 y^1 to be consistent with Code V.

prysm.polynomials.xy.xy_polynomial_sequence(mns, x, y, cartesian_grid=True)#

Contemporary XY monomial sequence.

Parameters:
  • mns (iterable of length 2 vectors) – sequence [(m1, n1), (m2, n2), …]

  • x (numpy.ndarray) – x coordinates

  • y (numpy.ndarray) – y coordinates

  • cartesian_grid (bool, optional) – if True, the input grid is assumed to be cartesian, i.e., x and y axes are aligned to the array dimensions arr[y,x] to accelerate the computation

Returns:

list of modes, in the same order as mns

Return type:

list

prysm.polynomials.xy.xy_polynomial(m, n, x, y, cartesian_grid=True)#

Contemporary XY monomial for a given m, n.

Parameters:
  • m (int) – x order

  • n (int) – y order

  • x (numpy.ndarray) – x coordinates

  • y (numpy.ndarray) – y coordinates

  • cartesian_grid (bool, optional) – if True, the input grid is assumed to be cartesian, i.e., x and y axes are aligned to the array dimensions arr[y,x] to accelerate the computation

Returns:

list of modes, in the same order as mns

Return type:

list

prysm.polynomials.xy.generalized_xy_polynomial_sequence(mns, x, y, seq_func, seq_func_kwargs=None, cartesian_grid=True)#

Generalized XY sequence.

Parameters:
  • mns (iterable of length 2 vectors) – sequence [(m1, n1), (m2, n2), …]

  • x (numpy.ndarray) – x coordinates

  • y (numpy.ndarray) – y coordinates

  • seq_func (callable) – signature seq_func(ns, x, **kwargs) -> list[numpy.ndarray] a function to generate a sequence of polynomials for one of the dimensions for example, cheby1_seq, legendre_seq, hermite_xx_seq, … from prysm.polynomials all work

  • seq_func_kwargs (dict, optional) – keyword arguments for seq_func

  • cartesian_grid (bool, optional) – if True, the input grid is assumed to be cartesian, i.e., x and y axes are aligned to the array dimensions arr[y,x] to accelerate the computation

Returns:

list of modes, in the same order as mns

Return type:

list

Q (Forbes)#

Tools for working with Q (Forbes) polynomials.

prysm.polynomials.qpoly.g_qbfs(n_minus_1)#

g(m-1) from oe-18-19-19700 eq. (A.15).

prysm.polynomials.qpoly.h_qbfs(n_minus_2)#

h(m-2) from oe-18-19-19700 eq. (A.14).

prysm.polynomials.qpoly.f_qbfs(n)#

f(m) from oe-18-19-19700 eq. (A.16).

prysm.polynomials.qpoly.Qbfs(n, x)#

Qbfs polynomial of order n at point(s) x.

Parameters:
  • n (int) – polynomial order

  • x (numpy.array) – point(s) at which to evaluate

Returns:

Qbfs_n(x)

Return type:

numpy.ndarray

prysm.polynomials.qpoly.change_basis_Qbfs_to_Pn(cs)#

Perform the change of basis from Qbfs to the auxiliary polynomial Pn.

The auxiliary polynomial is defined in A.4 of oe-18-19-19700 and is the shifted Chebyshev polynomials of the third kind.

Qbfs polynomials u^2(1-u^2)Qbfs_n(u^2) can be expressed as u^2(1-u^2)Pn(u^2) u in Forbes’ parlance is the normalized radial coordinate, so given points r in the range [0,1], use this function and then polynomials.cheby3(n, r*r). The u^2 (1 - u^2) is baked into the Qbfs function and will need to be applied by the caller for Cheby3.

Parameters:

cs (iterable) – sequence of polynomial coefficients, from order n=0..len(cs)-1

Returns:

array of same type as cs holding the coefficients that represent the same surface as a sum of shifted Chebyshev polynomials of the third kind

Return type:

numpy.ndarray

prysm.polynomials.qpoly.clenshaw_qbfs(cs, usq, alphas=None)#

Use Clenshaw’s method to compute a Qbfs surface from its coefficients.

Parameters:
  • cs (iterable of float) – coefficients for a Qbfs surface, from order 0..len(cs)-1

  • usq (numpy.ndarray) – radial coordinate(s) to evaluate, squared, notionally in the range [0,1] the variable u^2 from oe-18-19-19700

  • alphas (numpy.ndarray, optional) – array to store the alpha sums in, the surface is u^2(1-u^2) * (2 * (alphas[0]+alphas[1]) if not None, alphas should be of shape (len(s), x.shape) see _initialize_alphas if you desire more information

Returns:

Qbfs surface, the quantity u^2(1-u^2) S(u^2) from Eq. (3.13) note: excludes the division by phi, since c and rho are unknown

Return type:

numpy.ndarray

prysm.polynomials.qpoly.clenshaw_qbfs_der(cs, usq, j=1, alphas=None)#

Use Clenshaw’s method to compute Nth order derivatives of a sum of Qbfs polynomials.

Excludes base sphere and u^2(1-u^2) prefix

As an end-user, you are likely more interested in compute_zprime_Qbfs.

Parameters:
  • cs (iterable of float) – coefficients for a Qbfs surface, from order 0..len(cs)-1

  • usq (numpy.ndarray) – radial coordinate(s) to evaluate, squared, notionally in the range [0,1] the variable u^2 from oe-18-19-19700

  • j (int) – derivative order

  • alphas (numpy.ndarray, optional) –

    array to store the alpha sums in, if x = u * u, then S = (x * (1 - x)) * 2 * (alphas[0][0] + alphas[0][1]) S’ = … .. the same, but alphas[1][0] and alphas[1][1] S’’ = … … … … … … [2][0] … … ..[1][1] etc

    if not None, alphas should be of shape (j+1, len(cs), *x.shape) see _initialize_alphas if you desire more information

Returns:

the alphas array

Return type:

numpy.ndarray

prysm.polynomials.qpoly.product_rule(u, v, du, dv)#

The product rule of calculus, d/dx uv = u dv v du.

prysm.polynomials.qpoly.compute_z_zprime_Qbfs(coefs, u, usq)#

Compute the surface sag and first radial derivative of a Qbfs surface.

Excludes base sphere.

from Eq. 3.13 and 3.14 of oe-18-19-19700.

Parameters:
  • coefs (iterable) – surface coefficients for Q0..QN, N=len(coefs)-1

  • u (numpy.ndarray) – normalized radial coordinates (rho/rho_max)

  • usq (numpy.ndarray) – u^2

  • c (float) – best fit sphere curvature use c=0 for a flat base surface

Returns:

S, Sprime in Forbes’ parlance

Return type:

numpy.ndarray, numpy.ndarray

prysm.polynomials.qpoly.compute_z_zprime_Qcon(coefs, u, usq)#

Compute the surface sag and first radial derivative of a Qcon surface.

Excludes base sphere.

from Eq. 5.3 and 5.3 of oe-18-13-13851.

Parameters:
  • coefs (iterable) – surface coefficients for Q0..QN, N=len(coefs)-1

  • u (numpy.ndarray) – normalized radial coordinates (rho/rho_max)

  • usq (numpy.ndarray) – u^2

Returns:

S, Sprime in Forbes’ parlance

Return type:

numpy.ndarray, numpy.ndarray

prysm.polynomials.qpoly.Qbfs_sequence(ns, x)#

Qbfs polynomials of orders ns at point(s) x.

Parameters:
  • ns (Iterable of int) – polynomial orders

  • x (numpy.array) – point(s) at which to evaluate

Returns:

yielding one order of ns at a time

Return type:

generator of numpy.ndarray

prysm.polynomials.qpoly.Qcon(n, x)#

Qcon polynomial of order n at point(s) x.

Parameters:
  • n (int) – polynomial order

  • x (numpy.array) – point(s) at which to evaluate

Returns:

Qcon_n(x)

Return type:

numpy.ndarray

Notes

The argument x is notionally uniformly spaced 0..1. The Qcon polynomials are obtained by computing c = x^4. A transformation is then made, x => 2x^2 - 1 and the Qcon polynomials are defined as the jacobi polynomials with alpha=0, beta=4, the same order n, and the transformed x. The result of that is multiplied by c to yield a Qcon polynomial. Sums can more quickly be calculated by deferring the multiplication by c.

prysm.polynomials.qpoly.Qcon_sequence(ns, x)#

Qcon polynomials of orders ns at point(s) x.

Parameters:
  • ns (Iterable of int) – polynomial orders

  • x (numpy.array) – point(s) at which to evaluate

Returns:

yielding one order of ns at a time

Return type:

generator of numpy.ndarray

prysm.polynomials.qpoly.abc_q2d(n, m)#

A, B, C terms for 2D-Q polynomials. oe-20-3-2483 Eq. (A.3).

Parameters:
  • n (int) – radial order

  • m (int) – azimuthal order

Returns:

A, B, C

Return type:

float, float, float

prysm.polynomials.qpoly.G_q2d(n, m)#

G term for 2D-Q polynomials. oe-20-3-2483 Eq. (A.15).

Parameters:
  • n (int) – radial order

  • m (int) – azimuthal order

Returns:

G

Return type:

float

prysm.polynomials.qpoly.F_q2d(n, m)#

F term for 2D-Q polynomials. oe-20-3-2483 Eq. (A.13).

Parameters:
  • n (int) – radial order

  • m (int) – azimuthal order

Returns:

F

Return type:

float

prysm.polynomials.qpoly.g_q2d(n, m)#

Lowercase g term for 2D-Q polynomials. oe-20-3-2483 Eq. (A.18a).

Parameters:
  • n (int) – radial order less one (n - 1)

  • m (int) – azimuthal order

Returns:

g

Return type:

float

prysm.polynomials.qpoly.f_q2d(n, m)#

Lowercase f term for 2D-Q polynomials. oe-20-3-2483 Eq. (A.18b).

Parameters:
  • n (int) – radial order

  • m (int) – azimuthal order

Returns:

f

Return type:

float

prysm.polynomials.qpoly.Q2d(n, m, r, t)#

2D Q polynomial, aka the Forbes polynomials.

Parameters:
  • n (int) – radial polynomial order

  • m (int) – azimuthal polynomial order

  • r (numpy.ndarray) – radial coordinate, slope orthogonal in [0,1]

  • t (numpy.ndarray) – azimuthal coordinate, radians

Returns:

array containing Q2d_n^m(r,t) the leading coefficient u^m or u^2 (1 - u^2) and sines/cosines are included in the return

Return type:

numpy.ndarray

prysm.polynomials.qpoly.Q2d_sequence(nms, r, t)#

Sequence of 2D-Q polynomials.

Parameters:
  • nms (iterable of tuple) – (n,m) for each desired term

  • r (numpy.ndarray) – radial coordinates

  • t (numpy.ndarray) – azimuthal coordinates

Returns:

yields one term for each element of nms

Return type:

generator

prysm.polynomials.qpoly.change_of_basis_Q2d_to_Pnm(cns, m)#

Perform the change of basis from Q_n^m to the auxiliary polynomial P_n^m.

The auxiliary polynomial is defined in A.1 of oe-20-3-2483 and is the an unconventional variant of Jacobi polynomials.

For terms where m=0, see change_basis_Qbfs_to_Pn. This function only concerns those terms within the sum u^m a_n^m cos(mt) + b_n^m sin(mt) Q_n^m(u^2) sum

Parameters:
  • cns (iterable) – sequence of polynomial coefficients, from order n=0..len(cs)-1 and a given m (not |m|, but m, i.e. either “-2” or “+2” but not both)

  • m (int) – azimuthal order

Returns:

array of same type as cs holding the coefficients that represent the same surface as a sum of shifted Chebyshev polynomials of the third kind

Return type:

numpy.ndarray

prysm.polynomials.qpoly.abc_q2d_clenshaw(n, m)#

Special twist on A.3 for B.7.

prysm.polynomials.qpoly.clenshaw_q2d(cns, m, usq, alphas=None)#

Use Clenshaw’s method to compute the alpha sums for a piece of a Q2D surface.

Parameters:
  • cns (iterable of float) – coefficients for a Qbfs surface, from order 0..len(cs)-1

  • m (int) – azimuthal order for the cns

  • usq (numpy.ndarray) – radial coordinate(s) to evaluate, squared, notionally in the range [0,1] the variable u^2 from oe-18-19-19700

  • alphas (numpy.ndarray, optional) – array to store the alpha sums in, the surface is u^2(1-u^2) * (2 * (alphas[0]+alphas[1]) if not None, alphas should be of shape (len(s), *x.shape) see _initialize_alphas if you desire more information

Returns:

array containing components to compute the surface sag sum(cn Qn) = .5 alphas[0] - 2/5 alphas[3], if m=1 and N>2,

.5 alphas[0], otherwise

Return type:

alphas

prysm.polynomials.qpoly.clenshaw_q2d_der(cns, m, usq, j=1, alphas=None)#

Use Clenshaw’s method to compute Nth order derivatives of a Q2D surface.

This function is to be consumed by the other parts of prysm, and simply does the “alphas” computations (B.10) and adjacent Eqns

See compute_zprime_Q2D for this calculation integrated

Parameters:
  • cns (iterable of float) – coefficients for a Qbfs surface, from order 0..len(cs)-1

  • m (int) – azimuthal order

  • usq (numpy.ndarray) – radial coordinate(s) to evaluate, squared, notionally in the range [0,1] the variable u from oe-18-19-19700

  • j (int) – derivative order

  • alphas (numpy.ndarray, optional) – array to store the alpha sums in, if not None, alphas should be of shape (j+1, len(cs), *x.shape) see _initialize_alphas if you desire more information

Returns:

the alphas array

Return type:

numpy.ndarray

prysm.polynomials.qpoly.compute_z_zprime_Q2d(cm0, ams, bms, u, t)#

Compute the surface sag and first radial and azimuthal derivative of a Q2D surface.

Excludes base sphere.

from Eq. 2.2 and Appendix B of oe-20-3-2483.

Parameters:
  • cm0 (iterable) – surface coefficients when m=0 (inside curly brace, top line, Eq. B.1) span n=0 .. len(cms)-1 and mus tbe fully dense

  • ams (iterable of iterables) – ams[0] are the coefficients for the m=1 cosine terms, ams[1] for the m=2 cosines, and so on. Same order n rules as cm0

  • bms (iterable of iterables) –

    same as ams, but for the sine terms ams and bms must be the same length - that is, if an azimuthal order m is presnet in ams, it must be present in bms. The azimuthal orders need not have equal radial expansions.

    For example, if ams extends to m=3, then bms must reach m=3 but, if the ams for m=3 span n=0..5, it is OK for the bms to span n=0..3, or any other value, even just [0].

  • u (numpy.ndarray) – normalized radial coordinates (rho/rho_max)

  • t (numpy.ndarray) – azimuthal coordinate, in the range [0, 2pi]

Returns:

surface sag, radial derivative of sag, azimuthal derivative of sag

Return type:

numpy.ndarray, numpy.ndarray, numpy.ndarray

prysm.polynomials.qpoly.Q2d_nm_c_to_a_b(nms, coefs)#

Re-structure Q2D coefficients to the form needed by compute_z_zprime_Q2d.

Parameters:
  • nms (iterable) – sequence of [(n1, m1), (n2, m2), …] negative m encodes “sine term” while positive m encodes “cosine term”

  • coefs (iterable) – same length as nms, coefficients for mode n_m

Returns:

list 1 is cms, the “Qbfs” coefficients (m=0) list 2 is the “a” coefficients (cosine terms) list 3 is the “b” coefficients (sine terms)

lists 2 and 3 are lists-of-lists and begin from m=1 to m=M, containing an empty list if that order was not present in the input

Return type:

list, list, list

Jacobi#

High performance / recursive jacobi polynomial calculation.

prysm.polynomials.jacobi.weight(alpha, beta, x)#

The weight function of the jacobi polynomials for a given alpha, beta value.

prysm.polynomials.jacobi.recurrence_abc(n, alpha, beta)#

See A&S online - https://dlmf.nist.gov/18.9 .

Pn = (an-1 x + bn-1) Pn-1 - cn-1 * Pn-2

This function makes a, b, c for the given n, i.e. to get a(n-1), do recurrence_abc(n-1)

prysm.polynomials.jacobi.jacobi(n, alpha, beta, x)#

Jacobi polynomial of order n with weight parameters alpha and beta.

Parameters:
  • n (int) – polynomial order

  • alpha (float) – first weight parameter

  • beta (float) – second weight parameter

  • x (numpy.ndarray) – x coordinates to evaluate at

Returns:

jacobi polynomial evaluated at the given points

Return type:

numpy.ndarray

prysm.polynomials.jacobi.jacobi_sequence(ns, alpha, beta, x)#

Jacobi polynomials of orders ns with weight parameters alpha and beta.

Parameters:
  • ns (iterable) – sorted polynomial orders to return, e.g. [1, 3, 5, 7, …]

  • alpha (float) – first weight parameter

  • beta (float) – second weight parameter

  • x (numpy.ndarray) – x coordinates to evaluate at

Returns:

equivalent to array of shape (len(ns), len(x))

Return type:

generator

prysm.polynomials.jacobi.jacobi_der(n, alpha, beta, x)#

First derivative of Pn with respect to x, at points x.

Parameters:
  • n (int) – polynomial order

  • alpha (float) – first weight parameter

  • beta (float) – second weight parameter

  • x (numpy.ndarray) – x coordinates to evaluate at

Returns:

jacobi polynomial evaluated at the given points

Return type:

numpy.ndarray

prysm.polynomials.jacobi.jacobi_der_sequence(ns, alpha, beta, x)#

First partial derivative of Pn w.r.t. x for order ns, i.e. P_n’.

Parameters:
  • ns (iterable) – sorted orders to return, e.g. [1, 2, 3, 10] returns P1’, P2’, P3’, P10’

  • alpha (float) – first weight parameter

  • beta (float) – second weight parameter

  • x (numpy.ndarray) – x coordinates to evaluate at

Returns:

equivalent to array of shape (len(ns), len(x))

Return type:

generator

prysm.polynomials.jacobi.jacobi_sum_clenshaw(s, alpha, beta, x, alphas=None)#

Compute a weighted sum of Jacobi polynomials using Clenshaw’s method.

Parameters:
  • s (iterable) – weights in ascending order, beginning with P0, then P1, etc. must be fully dense when iterated

  • alpha (float) – first Jacobi shape parameter

  • beta (float) – second Jacobi shape parameter

  • x (numpy.ndarray or float_like) – coordinates to evaluate the sum at, orthogonal over [-1,1]

  • alphas (numpy.ndarray, optional) – array to store the alpha sums in, alphas[0] contains the sum and is returned if not None, alphas should be of shape (len(s), x.shape) see _initialize_alphas if you desire more information

Returns:

weighted sum of Jacobi polynomials

Return type:

numpy.ndarray

prysm.polynomials.jacobi.jacobi_sum_clenshaw_der(s, alpha, beta, x, j=1, alphas=None)#

Compute a weighted sum of partial derivatives w.r.t. x of Jacobi polynomials using Clenshaw’s method.

Notes

If the polynomial values and their derivatives are desired, pass alphas instead of leaving it None. alphas[0,0] will contain the sum of the polynomials, alphas[1,0] the sum of the first derivative, and so on.

Parameters:
  • s (iterable) – weights in ascending order, beginning with P0, then P1, etc. must be fully dense when iterated

  • alpha (float) – first Jacobi shape parameter

  • beta (float) – second Jacobi shape parameter

  • x (numpy.ndarray or float_like) – coordinates to evaluate the sum at, orthogonal over [-1,1]

  • j (int) – derivative order to compute

  • alphas (numpy.ndarray, optional) –

    array to store the alpha sums in, alphas[n] is the nth order derivative alpha terms with n=0 being the non-derivative terms.

    for a given n, the value of alphas[0] is the nth derivative of the surface sum if not None, alphas should be of shape (j+1, len(s), *x.shape) see _initialize_alphas if you desire more information

Returns:

alphas array, see alphas parameter documentation for meaning

Return type:

numpy.ndarray

Chebyshev#

Chebyshev polynomials.

prysm.polynomials.cheby.cheby1(n, x)#

Chebyshev polynomial of the first kind of order n.

Parameters:
  • n (int) – order to evaluate

  • x (numpy.ndarray) – point(s) at which to evaluate, orthogonal over [-1,1]

prysm.polynomials.cheby.cheby1_sequence(ns, x)#

Chebyshev polynomials of the first kind of orders ns.

Faster than chevy1 in a loop.

Parameters:
  • ns (Iterable of int) – orders to evaluate

  • x (numpy.ndarray) – point(s) at which to evaluate, orthogonal over [-1,1]

prysm.polynomials.cheby.cheby1_der(n, x)#

Partial derivative w.r.t. x of Chebyshev polynomial of the first kind of order n.

Parameters:
  • n (int) – order to evaluate

  • x (numpy.ndarray) – point(s) at which to evaluate, orthogonal over [-1,1]

prysm.polynomials.cheby.cheby1_der_sequence(ns, x)#

Partial derivative w.r.t. x of Chebyshev polynomials of the first kind of orders ns.

Faster than chevy1_der in a loop.

Parameters:
  • ns (Iterable of int) – orders to evaluate

  • x (numpy.ndarray) – point(s) at which to evaluate, orthogonal over [-1,1]

prysm.polynomials.cheby.cheby2(n, x)#

Chebyshev polynomial of the second kind of order n.

Parameters:
  • n (int) – order to evaluate

  • x (numpy.ndarray) – point(s) at which to evaluate, orthogonal over [-1,1]

prysm.polynomials.cheby.cheby2_sequence(ns, x)#

Chebyshev polynomials of the second kind of orders ns.

Faster than chevy1 in a loop.

Parameters:
  • ns (Iterable of int) – orders to evaluate

  • x (numpy.ndarray) – point(s) at which to evaluate, orthogonal over [-1,1]

prysm.polynomials.cheby.cheby2_der(n, x)#

Partial derivative w.r.t. x of Chebyshev polynomial of the second kind of order n.

Parameters:
  • n (int) – order to evaluate

  • x (numpy.ndarray) – point(s) at which to evaluate, orthogonal over [-1,1]

prysm.polynomials.cheby.cheby2_der_sequence(ns, x)#

Partial derivative w.r.t. x of Chebyshev polynomials of the second kind of orders ns.

Faster than chevy2_der in a loop.

Parameters:
  • ns (Iterable of int) – orders to evaluate

  • x (numpy.ndarray) – point(s) at which to evaluate, orthogonal over [-1,1]

prysm.polynomials.cheby.cheby3(n, x)#

Chebyshev polynomial of the third kind of order n.

Parameters:
  • n (int) – order to evaluate

  • x (numpy.ndarray) – point(s) at which to evaluate, orthogonal over [-1,1]

prysm.polynomials.cheby.cheby3_sequence(ns, x)#

Chebyshev polynomials of the third kind of orders ns.

Faster than chevy1 in a loop.

Parameters:
  • ns (Iterable of int) – orders to evaluate

  • x (numpy.ndarray) – point(s) at which to evaluate, orthogonal over [-1,1]

prysm.polynomials.cheby.cheby3_der(n, x)#

Partial derivative w.r.t. x of Chebyshev polynomial of the third kind of order n.

Parameters:
  • n (int) – order to evaluate

  • x (numpy.ndarray) – point(s) at which to evaluate, orthogonal over [-1,1]

prysm.polynomials.cheby.cheby3_der_sequence(ns, x)#

Partial derivative w.r.t. x of Chebyshev polynomials of the third kind of orders ns.

Faster than chevy1_der in a loop.

Parameters:
  • ns (Iterable of int) – orders to evaluate

  • x (numpy.ndarray) – point(s) at which to evaluate, orthogonal over [-1,1]

prysm.polynomials.cheby.cheby4(n, x)#

Chebyshev polynomial of the fourth kind of order n.

Parameters:
  • n (int) – order to evaluate

  • x (numpy.ndarray) – point(s) at which to evaluate, orthogonal over [-1,1]

prysm.polynomials.cheby.cheby4_sequence(ns, x)#

Chebyshev polynomials of the fourth kind of orders ns.

Faster than chevy1 in a loop.

Parameters:
  • ns (Iterable of int) – orders to evaluate

  • x (numpy.ndarray) – point(s) at which to evaluate, orthogonal over [-1,1]

prysm.polynomials.cheby.cheby4_der(n, x)#

Partial derivative w.r.t. x of Chebyshev polynomial of the fourth kind of order n.

Parameters:
  • n (int) – order to evaluate

  • x (numpy.ndarray) – point(s) at which to evaluate, orthogonal over [-1,1]

prysm.polynomials.cheby.cheby4_der_sequence(ns, x)#

Partial derivative w.r.t. x of Chebyshev polynomials of the fourth kind of orders ns.

Faster than chevy1_der in a loop.

Parameters:
  • ns (Iterable of int) – orders to evaluate

  • x (numpy.ndarray) – point(s) at which to evaluate, orthogonal over [-1,1]

Legendre#

Legendre polynomials.

prysm.polynomials.legendre.legendre(n, x)#

Legendre polynomial of order n.

Parameters:
  • n (int) – order to evaluate

  • x (numpy.ndarray) – point(s) at which to evaluate, orthogonal over [-1,1]

prysm.polynomials.legendre.legendre_sequence(ns, x)#

Legendre polynomials of orders ns.

Faster than legendre in a loop.

Parameters:
  • ns (int) – orders to evaluate

  • x (numpy.ndarray) – point(s) at which to evaluate, orthogonal over [-1,1]

prysm.polynomials.legendre.legendre_der(n, x)#

Partial derivative w.r.t. x of Legendre polynomial of order n.

Parameters:
  • n (int) – order to evaluate

  • x (numpy.ndarray) – point(s) at which to evaluate, orthogonal over [-1,1]

prysm.polynomials.legendre.legendre_der_sequence(ns, x)#

Partial derivative w.r.t. x of Legendre polynomials of orders ns.

Faster than legendre_der in a loop.

Parameters:
  • ns (int) – orders to evaluate

  • x (numpy.ndarray) – point(s) at which to evaluate, orthogonal over [-1,1]

Hermite#

Hermite Polynomials.

prysm.polynomials.hermite.hermite_He(n, x)#

Probabilist’s Hermite polynomial He of order n at points x.

Parameters:
  • n (int) – polynomial order

  • x (numpy.ndarray) – point(s) to evaluate at. Scalars and arrays both work.

Returns:

He_n(x)

Return type:

numpy.ndarray

prysm.polynomials.hermite.hermite_He_sequence(ns, x)#

Probabilist’s Hermite polynomials He of orders ns at points x.

Parameters:
  • ns (iterable of int) – rising polynomial orders, assumed to be sorted

  • x (numpy.ndarray) – point(s) to evaluate at. Scalars and arrays both work.

Returns:

equivalent to array of shape (len(ns), len(x))

Return type:

generator of numpy.ndarray

prysm.polynomials.hermite.hermite_He_der(n, x)#

First derivative of He_n with respect to x, at points x.

Parameters:
  • n (int) – polynomial order

  • x (numpy.ndarray) – point(s) to evaluate at. Scalars and arrays both work.

Returns:

d/dx[He_n(x)]

Return type:

numpy.ndarray

prysm.polynomials.hermite.hermite_He_der_sequence(ns, x)#

First derivative of He_[ns] with respect to x, at points x.

Parameters:
  • ns (iterable of int) – rising polynomial orders, assumed to be sorted

  • x (numpy.ndarray) – point(s) to evaluate at. Scalars and arrays both work.

Returns:

equivalent to array of shape (len(ns), len(x))

Return type:

generator of numpy.ndarray

prysm.polynomials.hermite.hermite_H(n, x)#

Physicist’s Hermite polynomial H of order n at points x.

Parameters:
  • n (int) – polynomial order

  • x (numpy.ndarray) – point(s) to evaluate at. Scalars and arrays both work.

Returns:

H_n(x)

Return type:

numpy.ndarray

prysm.polynomials.hermite.hermite_H_sequence(ns, x)#

Physicist’s Hermite polynomials H of orders ns at points x.

Parameters:
  • ns (iterable of int) – rising polynomial orders, assumed to be sorted

  • x (numpy.ndarray) – point(s) to evaluate at. Scalars and arrays both work.

Returns:

equivalent to array of shape (len(ns), len(x))

Return type:

generator of numpy.ndarray

prysm.polynomials.hermite.hermite_H_der(n, x)#

First derivative of H_n with respect to x, at points x.

Parameters:
  • n (int) – polynomial order

  • x (numpy.ndarray) – point(s) to evaluate at. Scalars and arrays both work.

Returns:

d/dx[H_n(x)]

Return type:

numpy.ndarray

prysm.polynomials.hermite.hermite_H_der_sequence(ns, x)#

First derivative of He_[ns] with respect to x, at points x.

Parameters:
  • ns (iterable of int) – rising polynomial orders, assumed to be sorted

  • x (numpy.ndarray) – point(s) to evaluate at. Scalars and arrays both work.

Returns:

equivalent to array of shape (len(ns), len(x))

Return type:

generator of numpy.ndarray

Dicksons#

Dickson Polynomials.

prysm.polynomials.dickson.dickson1(n, alpha, x)#

Dickson Polynomial of the first kind of order n.

Parameters:
  • n (int) – polynomial order

  • alpha (float) – shape parameter if alpha = -1, the dickson polynomials are Fibonacci Polynomials if alpha = 0, the dickson polynomials are the monomials x^n if alpha = 1, the dickson polynomials and cheby1 polynomials are related by D_n(2x) = 2T_n(x)

  • x (numpy.ndarray) – coordinates to evaluate the polynomial at

Returns:

D_n(x)

Return type:

numpy.ndarray

prysm.polynomials.dickson.dickson2(n, alpha, x)#

Dickson Polynomial of the second kind of order n.

Parameters:
  • n (int) – polynomial order

  • alpha (float) – shape parameter if alpha = -1, the dickson polynomials are Lucas Polynomials

  • x (numpy.ndarray) – coordinates to evaluate the polynomial at

Returns:

E_n(x)

Return type:

numpy.ndarray

prysm.polynomials.dickson.dickson1_sequence(ns, alpha, x)#

Sequence of Dickson Polynomial of the first kind of orders ns.

Parameters:
  • ns (iterable of int) – rising polynomial orders, assumed to be sorted

  • alpha (float) – shape parameter if alpha = -1, the dickson polynomials are Fibonacci Polynomials if alpha = 0, the dickson polynomials are the monomials x^n if alpha = 1, the dickson polynomials and cheby1 polynomials are related by D_n(2x) = 2T_n(x)

  • x (numpy.ndarray) – coordinates to evaluate the polynomial at

Returns:

equivalent to array of shape (len(ns), len(x))

Return type:

generator of numpy.ndarray

prysm.polynomials.dickson.dickson2_sequence(ns, alpha, x)#

Sequence of Dickson Polynomial of the second kind of orders ns.

Parameters:
  • ns (iterable of int) – rising polynomial orders, assumed to be sorted

  • alpha (float) – shape parameter if alpha = -1, the dickson polynomials are Lucas Polynomials

  • x (numpy.ndarray) – coordinates to evaluate the polynomial at

Returns:

D_n(x)

Return type:

numpy.ndarray