sgpykit.tools.polynomials_functions
Functions
|
Evaluate the k-th Chebyshev polynomial of the first kind on [a,b]. |
|
Evaluate the multidimensional Chebyshev polynomial of the first kind of order k (multi-index) on [a,b]^N. |
|
Evaluate the k-th generalized Laguerre polynomial orthonormal in [0, +inf). |
|
Evaluate the multidimensional generalized Laguerre polynomial of order k (multi-index) orthonormal on [0, +inf)^N with respect to the weight function: rho = prod_i beta_i^(alpha_i+1)/Gamma(alpha_i+1) * x^alpha_i * exp(-beta_i * x), where alpha_i > -1 and beta_i > 0. |
|
Evaluate the k-th orthonormal Hermite polynomial at points x. |
|
Evaluate the multidimensional Hermite polynomial of order k (multi-index) orthonormal on [-inf, +inf]^N with respect to the Gaussian weight function. |
|
Evaluate the k-th Jacobi probabilistic polynomial orthonormal in [a,b]. |
|
Evaluate the multidimensional probabilistic Jacobi polynomial of order k (multi-index) orthonormal on [a1,b1] x [a2,b2] x [a3,b3] x . |
|
Evaluate the Lagrange basis polynomial at non-grid points. |
|
Evaluate the Lagrange basis polynomial at non-grid points. |
|
Evaluate the multidimensional Lagrange function centered at current_knot in non_grid_points. |
|
Evaluate the k-th Laguerre polynomial orthonormal in [0, +inf) w.r.t. |
|
Evaluate the multidimensional Laguerre polynomial of order k (multi-index) orthonormal on [0, +inf)^N with respect to rho=prod_i lambda_i*exp(-lambda_i*x), lambda_i>0 on the list of points X (each column is a point in R^N). |
|
Evaluate the k-th orthonormal Legendre polynomial on the interval [a, b]. |
|
Evaluate the multidimensional Legendre polynomial of order k (multi-index) orthonormal on [a,b]^N on the list of points X (each column is a point in R^N). |
|
Evaluate the k-th standard Legendre polynomial. |
|
Interpolates values at non_grid_points using a Lagrange interpolant built from given grid_points and tabulated function_on_grid. |
- sgpykit.tools.polynomials_functions.cheb_eval(x, k, a, b)[source]
Evaluate the k-th Chebyshev polynomial of the first kind on [a,b].
- Parameters:
- xarray_like
Points at which to evaluate the Chebyshev polynomial.
- kint
Degree of the Chebyshev polynomial.
- afloat
Lower bound of the interval.
- bfloat
Upper bound of the interval.
- Returns:
- Lndarray
Values of the k-th Chebyshev polynomial evaluated at x.
- sgpykit.tools.polynomials_functions.cheb_eval_multidim(X, k, a, b)[source]
Evaluate the multidimensional Chebyshev polynomial of the first kind of order k (multi-index) on [a,b]^N.
Evaluates the multidim Chebyshev polynomial of the first kind of order k (multi-index) on [a,b]^N on the list of points X (each column is a point in R^N). a,b are scalar values a,b may differ in each direction, so that actually the domain is not [a,b]^N, but [a1,b1] x [a2,b2] x [a3,b3] x … [aN,bN] in this case, a,b are defined as vectors, each with N components
- Parameters:
- Xnumpy.ndarray
Array of points, shape (N, nb_pts), where each column is a point in R^N.
- knumpy.ndarray
Multi-index of the polynomial order, shape (N,).
- anumpy.ndarray or float
Lower bounds of the domain. If scalar, it is expanded to a vector of length N.
- bnumpy.ndarray or float
Upper bounds of the domain. If scalar, it is expanded to a vector of length N.
- Returns:
- Cnumpy.ndarray
Evaluated polynomial values at each point, shape (nb_pts,).
- sgpykit.tools.polynomials_functions.generalized_lagu_eval(x, k, alpha, beta)[source]
Evaluate the k-th generalized Laguerre polynomial orthonormal in [0, +inf).
This function returns the values of the k-th generalized Laguerre polynomial orthonormal in [0, +inf) with respect to the weight function rho = beta^(alpha+1)/Gamma(alpha+1) * x^alpha * exp(-beta*x) at the points x. The polynomials start from k=0: L_0(x) = 1, L_1(x) = -x*beta + alpha + 1. This function expresses L as a function of the standard generalized Laguerre “probabilistic” polynomial (i.e., orthogonal w.r.t. rho = x^alpha * exp(-x)/Gamma(alpha+1)), which are recursively calculated through the function standard_generalized_lagu_eval.
- Parameters:
- xarray_like
Points at which to evaluate the polynomial. Can be a matrix.
- kint
Degree of the polynomial.
- alphafloat
Parameter of the generalized Laguerre polynomial.
- betafloat
Parameter of the generalized Laguerre polynomial.
- Returns:
- Lndarray
Values of the k-th generalized Laguerre polynomial at the points x.
- sgpykit.tools.polynomials_functions.generalized_lagu_eval_multidim(X, k, alpha, beta)[source]
Evaluate the multidimensional generalized Laguerre polynomial of order k (multi-index) orthonormal on [0, +inf)^N with respect to the weight function: rho = prod_i beta_i^(alpha_i+1)/Gamma(alpha_i+1) * x^alpha_i * exp(-beta_i * x), where alpha_i > -1 and beta_i > 0.
- Parameters:
- Xnumpy.ndarray
Array of points where the polynomial is evaluated. Each column is a point in R^N.
- knumpy.ndarray
Multi-index defining the polynomial order.
- alphanumpy.ndarray or scalar
Parameter alpha for the generalized Laguerre polynomial. Can be scalar or array.
- betanumpy.ndarray or scalar
Parameter beta for the generalized Laguerre polynomial. Can be scalar or array.
- Returns:
- Lnumpy.ndarray
Array of polynomial evaluations at the given points.
- sgpykit.tools.polynomials_functions.herm_eval(x, k, mi, sigma)[source]
Evaluate the k-th orthonormal Hermite polynomial at points x.
This function returns the values of the k-th Hermite polynomial orthonormal in (-inf, +inf) with respect to the weight function rho = 1/sqrt(2 pi sigma^2) * exp(-(x-mi)^2/(2*sigma^2)) at the points x (x can be a matrix as well).
The polynomials start from k=0: H_0(x) = 1, H_1(x) = (x - mi)/sigma. This function expresses H as a function of the standard Hermite “probabilistic” polynomial (i.e., orthogonal with respect to rho=1/sqrt(2 pi) * exp(-x^2/2)), which are recursively calculated through the function standard_herm_eval.
- Parameters:
- xarray_like
Points at which to evaluate the polynomial.
- kint
Degree of the Hermite polynomial.
- mifloat
Mean of the Gaussian distribution.
- sigmafloat
Standard deviation of the Gaussian distribution.
- Returns:
- Hndarray
Values of the k-th orthonormal Hermite polynomial at points x.
- sgpykit.tools.polynomials_functions.herm_eval_multidim(X, k, mu, sig)[source]
Evaluate the multidimensional Hermite polynomial of order k (multi-index) orthonormal on [-inf, +inf]^N with respect to the Gaussian weight function.
- Parameters:
- Xnumpy.ndarray
Array of points where the polynomial is evaluated. Each column is a point in R^N.
- knumpy.ndarray
Multi-index defining the order of the polynomial in each dimension.
- mufloat or numpy.ndarray
Mean of the Gaussian weight function. Can be a scalar or a vector of length N.
- sigfloat or numpy.ndarray
Standard deviation of the Gaussian weight function. Can be a scalar or a vector of length N.
- Returns:
- Hnumpy.ndarray
Array of polynomial evaluations at the points in X.
- sgpykit.tools.polynomials_functions.jacobi_prob_eval(x, k, alpha, beta, a, b)[source]
Evaluate the k-th Jacobi probabilistic polynomial orthonormal in [a,b].
This function returns the values of the k-th Jacobi “probabilistic” polynomial orthonormal in [a,b] with respect to the weight function: rho(x) = Gamma(alpha+beta+2)/(Gamma(alpha+1)*Gamma(beta+1)*(b-a)^(alpha+beta+1)) * (x-a)^alpha * (b-x)^beta in the points x (x can be a matrix as well).
Note that these are NOT the classical Jacobi polynomials (which are orthonormal in [a,b] w.r.t. (b-x)^alpha*(x-a)^beta).
- Parameters:
- xarray_like
Points at which to evaluate the polynomial.
- kint
Degree of the polynomial.
- alphafloat
Parameter of the Jacobi polynomial.
- betafloat
Parameter of the Jacobi polynomial.
- afloat
Lower bound of the interval.
- bfloat
Upper bound of the interval.
- Returns:
- Lndarray
Values of the k-th Jacobi probabilistic polynomial at points x.
Notes
The polynomials start from k=0: L_0(x) = 1, L_1(x) = 0.5*(alpha-beta) + 0.5*(alpha+beta+2)*(-2*x+a+b)/(b-a)
This function expresses L as a function of the standard Jacobi polynomial (i.e. orthoGONAL w.r.t. rho(x)=(1-x)^alpha*(1+x)^beta, which are recursively calculated through the function standard_jac_eval.
- sgpykit.tools.polynomials_functions.jacobi_prob_eval_multidim(X, k, alpha, beta, a, b)[source]
Evaluate the multidimensional probabilistic Jacobi polynomial of order k (multi-index) orthonormal on [a1,b1] x [a2,b2] x [a3,b3] x … [aN,bN] with respect to rho=prod_i Gamma(alpha_i+beta_i+2)/(Gamma(alpha_i+1)*Gamma(beta_i+1)*(b-a)^(alpha_i+beta_i+1))*(x-a)^alpha_i*(b-x)^beta_i, alpha_i,beta_i>-1 on the list of points X (each column is a point in R^N).
- Parameters:
- Xnumpy.ndarray
Array of points, shape (N, nb_pts), where N is the dimension and nb_pts is the number of points.
- knumpy.ndarray
Multi-index of the polynomial, shape (N,).
- alphanumpy.ndarray or float
Parameter alpha for the Jacobi polynomial, can be scalar or array of shape (N,).
- betanumpy.ndarray or float
Parameter beta for the Jacobi polynomial, can be scalar or array of shape (N,).
- anumpy.ndarray or float
Lower bound of the interval, can be scalar or array of shape (N,).
- bnumpy.ndarray or float
Upper bound of the interval, can be scalar or array of shape (N,).
- Returns:
- Lnumpy.ndarray
Evaluated polynomial values, shape (nb_pts,).
- sgpykit.tools.polynomials_functions.lagr_eval(current_knot, other_knots, non_grid_points)[source]
Evaluate the Lagrange basis polynomial at non-grid points.
Builds the Lagrange function L(x) such that: L(current_knot) = 1 L(other_knots) = 0
and returns L = L(non_grid_points)
- Parameters:
- current_knotfloat
The knot where the Lagrange function should be 1.
- other_knotsarray_like
A list or array of knots where the Lagrange function should be 0.
- non_grid_pointsarray_like
The points at which to evaluate the Lagrange function.
- Returns:
- Lndarray
The evaluated Lagrange function at non_grid_points.
- sgpykit.tools.polynomials_functions.lagr_eval_fast(current_knot, other_knots, ok_len, non_grid_points, ng_size)[source]
Evaluate the Lagrange basis polynomial at non-grid points.
This function builds the Lagrange function L(x) such that L(current_knot)=1 and L(other_knots)=0, and returns L evaluated at non_grid_points.
- Parameters:
- current_knotfloat
The knot where the Lagrange polynomial equals 1.
- other_knotsarray_like
The knots where the Lagrange polynomial equals 0.
- ok_lenint
The length of other_knots.
- non_grid_pointsarray_like
The points at which to evaluate the Lagrange polynomial.
- ng_sizeint or tuple
The size of non_grid_points. If scalar, it is treated as the size of a square matrix. If tuple, it is treated as the shape of non_grid_points.
- Returns:
- Lndarray
The evaluated Lagrange polynomial at non_grid_points.
- sgpykit.tools.polynomials_functions.lagr_eval_multidim(current_knot, knots_per_dim, non_grid_points)[source]
Evaluate the multidimensional Lagrange function centered at current_knot in non_grid_points.
- Parameters:
- current_knotarray_like
The knot where the Lagrange function is centered.
- knots_per_dimlist or array_like
List of arrays, where each array contains all the knots in each dimension.
- non_grid_pointsarray_like
Matrix where each row is a point to be evaluated.
- Returns:
- Lndarray
The evaluated multidimensional Lagrange function at non_grid_points.
- sgpykit.tools.polynomials_functions.lagu_eval(x, k, lambda_)[source]
Evaluate the k-th Laguerre polynomial orthonormal in [0, +inf) w.r.t. rho=lambda*exp(-lambda*x).
This function computes the values of the k-th Laguerre polynomial, which is orthonormal with respect to the weight function rho = lambda * exp(-lambda * x) on the interval [0, +inf). The polynomials are defined starting from k=0, with L_0(x) = 1 and L_1(x) = -x*lambda + 1.
- Parameters:
- xarray_like
Points at which to evaluate the polynomial (can be a matrix).
- kint
Degree of the Laguerre polynomial.
- lambda_float
Parameter of the weight function (must be positive).
- Returns:
- Lndarray
Values of the k-th Laguerre polynomial at the points x.
- sgpykit.tools.polynomials_functions.lagu_eval_multidim(X, k, lambda_)[source]
Evaluate the multidimensional Laguerre polynomial of order k (multi-index) orthonormal on [0, +inf)^N with respect to rho=prod_i lambda_i*exp(-lambda_i*x), lambda_i>0 on the list of points X (each column is a point in R^N). Lambda can be a scalar value.
- Parameters:
- Xnumpy.ndarray
Array of points, shape (N, nb_pts), where N is the dimension and nb_pts is the number of points.
- knumpy.ndarray
Multi-index of the polynomial, shape (N,).
- lambda_float or numpy.ndarray
Scale parameter(s) for the Laguerre polynomial. If scalar, it is expanded to a vector of length N.
- Returns:
- Lnumpy.ndarray
Evaluated polynomial values, shape (nb_pts,).
- sgpykit.tools.polynomials_functions.lege_eval(x, k, a, b)[source]
Evaluate the k-th orthonormal Legendre polynomial on the interval [a, b].
This function computes the values of the k-th Legendre polynomial, which is orthonormal with respect to the uniform weight function rho = 1/(b-a), at the points specified by x. The polynomials are defined starting from k=0: L_0(x) = 1, L_1(x) = (2*x - a - b)/(b - a).
The computation is performed by transforming the input points from the interval [a, b] to the standard interval [-1, 1], evaluating the standard Legendre polynomial (orthogonal with respect to rho=1) at these transformed points, and then applying the necessary normalizations to account for the interval [a, b].
- Parameters:
- xarray_like
The points at which to evaluate the polynomial. Can be a scalar, list, or ndarray.
- kint
The degree of the Legendre polynomial.
- afloat
The lower bound of the interval.
- bfloat
The upper bound of the interval.
- Returns:
- Lndarray
The values of the k-th orthonormal Legendre polynomial evaluated at x.
- sgpykit.tools.polynomials_functions.lege_eval_multidim(X, k, a, b)[source]
Evaluate the multidimensional Legendre polynomial of order k (multi-index) orthonormal on [a,b]^N on the list of points X (each column is a point in R^N). a,b are scalar values.
- Parameters:
- Xnumpy.ndarray
Array of points, shape (N, nb_pts), where N is the dimension and nb_pts is the number of points.
- knumpy.ndarray
Multi-index of the polynomial, shape (N,).
- aint or numpy.ndarray
Lower bound(s) of the domain. If scalar, it is broadcasted to all dimensions.
- bint or numpy.ndarray
Upper bound(s) of the domain. If scalar, it is broadcasted to all dimensions.
- Returns:
- Lnumpy.ndarray
Evaluated polynomial values, shape (nb_pts,).
Notes
a,b may differ in each direction, so that actually the domain is not [a,b]^N, but [a1,b1] x [a2,b2] x [a3,b3] x … [aN,bN] in this case, a,b are defined as vectors, each with N components.
- sgpykit.tools.polynomials_functions.standard_lege_eval(x, k)[source]
Evaluate the k-th standard Legendre polynomial.
This function computes the values of the k-th standard Legendre polynomial, which is orthogonal (but not orthonormal) with respect to the uniform weight function rho=1 on the interval [-1, 1], at the points specified by x. The polynomials are defined starting from k=0: L_0(x) = 1, L_1(x) = x.
- Parameters:
- xarray_like
The points at which to evaluate the polynomial. Can be a scalar, list, or ndarray.
- kint
The degree of the Legendre polynomial.
- Returns:
- Lndarray
The values of the k-th standard Legendre polynomial evaluated at x.
- sgpykit.tools.polynomials_functions.univariate_interpolant(grid_points, function_on_grid, non_grid_points)[source]
Interpolates values at non_grid_points using a Lagrange interpolant built from given grid_points and tabulated function_on_grid.
- Parameters:
- grid_points(M,) array_like
Interpolation nodes.
- function_on_grid(V, M) array_like
Values of V different functions sampled at the M grid points.
- non_grid_points(N,) array_like
Points where we want to evaluate the interpolation.
- Returns:
- f_values(V, N) ndarray
Interpolated values at non_grid_points.
Notes
UNIVARIATE_INTERPOLANT interpolates a univariate function on grid points, i.e. evaluates the Lagrange polynomial approximation on a generic point of the parameter space.
- F_VALUES = UNIVARIATE_INTERPOLANT(GRID_POINTS,FUNCTION_ON_GRID,NON_GRID_POINTS) evaluates the
lagrangian interpolant of a vector function F: R -> R^V based on the points contained in the vector GRID_POINTS.
GRID POINTS is a row vector with the nodes of interpolation. FUNCTION_ON_GRID is a matrix containing the evaluation of F on the points GRID_POINTS.
Its dimensions are V x length(GRID_POINTS).
NON_GRID_POINTS is a row vector of points where one wants to evaluate the polynomial approximation. F_VALUES is a matrix containing the evaluation of the function F in each of the NON_GRID_POINTS.
Its dimensions are V X length(NON_GRID_POINTS)