PCE Base Class#

class tesuract.PCEBuilder(order=1, customM=None, mindex_type='total_order', coef=None, a=None, b=None, polytype='Legendre', normalized=False, store_phi=False, input_range=None, use_sklearn_poly_features=False)[source]#

Base class for building a multivariate polynomial basis.

This class creates a multi-variate polynomial object, aka as a polynomial chaos model. The expansion looks like

\[\sum_{i=1}^N c_i \Phi_i(\mathbf{x})\]

where \(N\) is the number of polynomial terms, and \(\Phi_i:\mathbf{x} \in \mathbb{R}^d \mapsto \mathbb{R}\) is the multivariate basis function which takes the form

\[\Phi_i = \prod_{j=1}^d L_{\alpha_j^{(i)}}(x_i),\]

where \(\alpha_j\) is an integer tuple of size \(d\) which represents the multiindex, and \(L:\mathbb{R} \mapsto \mathbb{R}\) are the one-dimensional Legendre polynomials. So for example, in two dimensions, we can have \(\Phi_i = xy\). The construction of the multiindex is done behind the scenes in the multindex module.

Below we list the parameters to construct the object followed by the class attributes and returns.

Parameters
orderint, default=1

Description of the order of the polynomials in the expansion. For total order, the order is the maximum polynomial order for each basis function per dimension.

customMnumpy.ndarray, default=None

An integer numpy aray of size (# of basis functions, # of dimensions). Each row represents the basis function and each column represents the order of the 1d polynomial, \(L_i\).

mindex_type{‘total_order’, ‘hyperbolic’}, default=’total_order’

Different types of multiindex generators. Total order produces basis vectors that have a maximum order as defined by order. Hyperbolic order is similar but generates fewer cross terms.

coefnumpy.ndarray, default=None

1d numpy coefficient array if defining polynomial. Must be the same length as the number of basis elements, i.e. length of multiindex array.

polytype{‘Legendre’}, default=’Legendre’

A string representing the type of polynomial. So far we only include Legendre polynomial construction defined on [-1,1]

normalizedbool, default=False

Whether or not to use normalized polynomials such that \(\int_{-1}^{1}L_i(x)dx = 1\) or \(\frac{2}{2n+1}\) if True or False, respecively, where \(n\) is the order of the polynomial (only for Legendre polynomials).

Attributes
dimint

Dimension of the polynomial, defined after construction.

mindexndarray of shape (nbasis, dim)

Integer array of the multiindex which describes the structure of the polynomial basis.

nPCTermsint

The number of basis terms of the multi-variate polynomial.

coefndarray of shape (nbasis,)

Coefficient array

Notes

As of now, the base class requires the input range to be on the default range \([-1,1]\) for all dimensions. We have included a useful preprocessing utility (DomainScaler) to transform the domain easily to the canonical target range from any input range. In the future this will be an option in the this class.

Examples

>>> from tesuract import PCEBuilder
>>> p = PCEBuilder(order=3,normalized=True)
>>> print(p.mindex)

Methods

compile:

pre-processing for fit, e.g., multiindex generation

computeMoments:

compute means and variances using PCE coefficients

compile(dim)[source]#

Setup for instantiating the basis class

Constructs the multi-dimensional multi-index which defines the polynomial basis elements. Note that this is only done once during the fit method of the class, unless the mindex variable is undefined.

The multi-index array is of size \(N \times dim\) where \(N\) is the number of basis elements. The multi-index determine the order or degree of the univariate polynomial in

\[\Phi_i = \prod_{j=1}^d L_{\alpha_j^{(i)}}(x_i).\]

So for example, \(\alpha_j=2\) corresponds to a second order or quadratic Legendre polynomial.

Lastly, this method also defines the number of basis elements.

Parameters
dimint

dimension of the polynomials

Returns
self: object

Returns the object itself.

transform(X)[source]#

Transform the given \(X\) coordinates to the polynomial space.

This method essentially performs a high dimensional kernel mapping onto a polynomial space spanned by the Legendre polynomials. This is similar to sklearn’s PolynomialFeatures, except here the features are multi-variate Legendre; the big difference being that the features are uncorrelated.

Alternatively, the resulting

Parameters
Xnumpy.ndarray of shape (nsamples, dim)

feature matrix where each row is a sample of the feature space.

Returns
Phinumpy.ndarray of shape (nsamples, nPCTerms)

returns the polynomial feature map that transforms each sample in \(X\) of dimension \(d\) to dimension \(nPCTerms\).

fit_transform(X)[source]#

Fit and transform the given \(X\) coordinates to the polynomial space.

This method essentially performs a high dimensional kernel mapping onto a polynomial space spanned by the Legendre polynomials. This is similar to sklearn’s PolynomialFeatures, except here the features are multi-variate Legendre; the big difference being that the features are uncorrelated.

Alternatively, the resulting

Parameters
Xnumpy.ndarray of shape (nsamples, dim)

feature matrix where each row is a sample of the feature space.

Returns
Phinumpy.ndarray of shape (nsamples, nPCTerms)

returns the polynomial feature map that transforms each sample in \(X\) of dimension \(d\) to dimension \(nPCTerms\).

polyeval(X=None, c=None)[source]#

Method to evaluate the polynomial.

Evaluates the polynomial for a given set of coefficients and data points.

Parameters
Xnumpy.ndarray of shape (nsamples, dim)

feature matrix where each row is a sample of the feature space.

cnumpy.ndarray of shape (nPCTerms,) (optional)

coefficient array. The evaluation is simply np.dot(Phi,c). The coefficient array can also be internally set so that it does not need to be fed in each time we need to evaluate the polynomial.

Returns
yevalnumpy.ndarray of shape (X.shape[0],)

The scalar outputs of the multivariate polynomial evaluated at the feature matrix data points.

eval(X=None, c=None)[source]#

Duplicate of polyeval

computeSobol(c=None)[source]#

Compute Sobol total order variance based sensitivity indices

Depending on the Legendre polynomial type (normalized vs non-normalized) the formula will be different. The total order sensitivity is given by

\[S_i = \sum_{k} \gamma^{-2}_{\beta^i_k} c^2_{\beta^i_k}\]

where \(\{\beta^i_k\}_{k=\dots}\) are the indices that contain at least the \(i^{th}\) dimension, and \(\gamma\) is the square root norm of that particular basis polynomial (which is 1 for the normalized case).

For this method, we return the normalized Sobol indices, i.e.

\[T_i \doteq S_i/S_{T}\]

where \(S_{T}\) is the total variance, i.e. \(\{\beta^i_k\}_{k=\dots}\) are the entire set of basis functions. Thus, the total order sensitivity indices are always less than 1, although their sum can be greater than \(S_{T}\).

This is another advantage of using the Legendre polynomials, in that they are uncorrelated so that their feature importance calculations are very easy to compute.

Parameters
cnumpy.ndarray of shape (nPCTerms,)

coefficient to determine the Sobol indices.

Returns
Tnumpy.ndarray of shape (dim,)

The total order Sobol sensitivity indices for each dimension

computeMoments(c=None)[source]#

Methods to compute the mean and variance of the resulting polynomial

Assuming a uniform distribution over \([-1,1]\), this method computes the mean and variance. The mean is the coefficient of the \(0^{th}\) order term, while the variance is the weighted sum of squares of the other coefficients. The weighted sum is determined by using either the normalized or non-normalized Legendre polynomials. Finally, the weight is scaled by .5 correpsonding to the density function for the uniform distribution on the previously stated domain.

\[\mu = \frac{1}{2} c_0\]
\[\sigma^2 = \frac{1}{2} \sum_{i=1}^N \gamma^2_i c^2_i\]
Parameters
cnumpy.ndarray of shape (nPCTerms, )

coefficient of the polynomial basis. Must have the same number of elements as the number of basis elements.

Returns
mu, varfloat,float

mean and variance floating point values