Polynomial#

The polynomial parent class; one of the main building blocks in Effective Quadratures.

class equadratures.poly.Graphpolys(Graph, data_train, poly, edge_weight=15)[source]#

Constructor for generating polynomials over graphs.

Parameters
  • Graph (networkx.Graph) – A networkx Graph object.

  • data_train (dict) – A dictionary with training data.

  • poly (Poly) – A Poly object that defines a polynomial over the subspace.

  • edge_weight (float, optional) – Edge weight for Graph edges, default=15.

fit()[source]#

Fit polynomial to data.

predict(data, score=True)[source]#

Predict.

Parameters
  • data (dict) – A dictionary with test data.

  • score (bool, optional) – returns mean error along with prediction if true.

Returns

  • Y_pred (numpy.nparray) – Predicted data.

  • mean_error (numpy.ndarray) – Mean error in prediction.

class equadratures.poly.Poly(parameters, basis, method=None, sampling_args=None, solver_args={}, variable=None, override_cardinality=False)[source]#

Definition of a polynomial object.

Parameters
  • parameters (list) – A list of parameters, where each element of the list is an instance of the Parameter class.

  • basis (Basis) – An instance of the Basis class corresponding to the multi-index set used.

  • method (str, optional) – The method used for computing the coefficients. Should be one of: compressed-sensing, least-squares, minimum-norm, numerical-integration, least-squares-with-gradients, least-absolute-residual, huber, elastic-net, relevance-vector-machine or custom-solver. See Solver for further details.

  • sampling_args (dict, optional) –

    A dict containing optional arguments centered around the specific sampling strategy:

    • mesh (str): Avaliable options are: monte-carlo, sparse-grid, tensor-grid, induced, or user-defined. Note that when the sparse-grid option is invoked, the sparse pseudospectral approximation method [1] is the adopted. One can think of this as being the correct way to use sparse grids in the context of polynomial chaos [2] techniques.

    • subsampling-algorithm (str): The subsampling-algorithm input refers to the optimisation technique for subsampling. In the aforementioned four sampling strategies, we generate a logarithm factor of samples above the required amount and prune down the samples using an optimisation technique (see [1]). Existing optimisation strategies include: qr, lu, svd, newton. These refer to QR with column pivoting [2], LU with row pivoting [3], singular value decomposition with subset selection [2] and a convex relaxation via Newton’s method for determinant maximization [4]. Note that if the tensor-grid option is selected, then subsampling will depend on whether the Basis argument is a total order index set, hyperbolic basis or a tensor order index set.

    • sampling-ratio (float): Denotes the extent of undersampling or oversampling required. For values equal to unity (default), the number of rows and columns of the associated Vandermonde-type matrix are equal.

    • sample-points (numpy.ndarray): A numpy ndarray with shape (number_of_observations, dimensions) that corresponds to a set of sample points over the parameter space.

    • sample-outputs (numpy.ndarray): A numpy ndarray with shape (number_of_observations, 1) that corresponds to model evaluations at the sample points. Note that if sample-points is provided as an input, then the code expects sample-outputs too.

    • sample-gradients (numpy.ndarray): A numpy ndarray with shape (number_of_observations, dimensions) that corresponds to a set of sample gradient values over the parameter space.

  • solver_args (dict, optional) – Optional arguments centered around the specific solver used for computing the coefficients. See Solver.

  • override_cardinality (bool, optional) – Option to override the soft cardinality limit of 50E3. By default the limit is enforced in order to avoid excessively long run times when the polynomial coefficients are computed.

Examples

>>> # Subsampling from a tensor grid
>>> param = eq.Parameter(distribution='uniform', lower=-1., upper=1., order=3)
>>> basis = eq.Basis('total order')
>>> poly = eq.Poly(parameters=[param, param], basis=basis, method='least-squares' , sampling_args={'mesh':'tensor-grid', 'subsampling-algorithm':'svd', 'sampling-ratio':1.0})
>>> # User-defined data with compressive sensing
>>> X = np.loadtxt('inputs.txt')
>>> y = np.loadtxt('outputs.txt')
>>> param = eq.Parameter(distribution='uniform', lower=-1., upper=1., order=3)
>>> basis = eq.Basis('total order')
>>> poly = eq.Poly([param, param], basis, method='compressive-sensing', sampling_args={'sample-points':X_red,'sample-outputs':Y_red})
>>> # Using a sparse grid
>>> param = eq.Parameter(distribution='uniform', lower=-1., upper=1., order=3)
>>> basis = eq.Basis('sparse-grid', level=7, growth_rule='exponential')
>>> poly = eq.Poly(parameters=[param, param], basis=basis, method='numerical-integration')

References

  1. Constantine, P. G., Eldred, M. S., Phipps, E. T., (2012) Sparse Pseudospectral Approximation Method. Computer Methods in Applied Mechanics and Engineering. 1-12. Paper

  2. Xiu, D., Karniadakis, G. E., (2002) The Wiener-Askey Polynomial Chaos for Stochastic Differential Equations. SIAM Journal on Scientific Computing, 24(2), Paper

  3. Seshadri, P., Iaccarino, G., Ghisu, T., (2018) Quadrature Strategies for Constructing Polynomial Approximations. Uncertainty Modeling for Engineering Applications. Springer, Cham, 2019. 1-25. Preprint

  4. Seshadri, P., Narayan, A., Sankaran M., (2017) Effectively Subsampled Quadratures for Least Squares Polynomial Approximations. SIAM/ASA Journal on Uncertainty Quantification, 5(1). Paper

  5. Rogers, S., Girolami, M., (2016) Variability in predictions. In: A First Course in Machine Learning, Second Edition (2nd. ed.). Chapman & Hall/CRC. Book

get_coefficients()[source]#

Returns the coefficients of the polynomial approximation.

Returns

The coefficients with size (number_of_coefficients, 1).

Return type

numpy.ndarray

get_conditional_kurtosis_indices(order)[source]#

Computes the kurtosis indices.

Parameters

order (int) – The highest order of the kurtosis indices required.

Returns

A dict comprising of kurtosis indices and constitutent mixed orders of the parameters.

Return type

dict

get_conditional_skewness_indices(order)[source]#

Computes the skewness indices.

Parameters

order (int) – The highest order of the skewness indices required.

Returns

A dict comprising of skewness indices and constitutent mixed orders of the parameters.

Return type

dict

get_mean_and_variance()[source]#

Computes the mean and variance of the model.

Returns

Tuple containing two floats; the approximated mean and variance from the polynomial fit.

Return type

tuple

get_model_evaluations()[source]#

Returns the output model evaluations used to fit the polynomial.

Returns

Array containing the output model evaluations.

Return type

numpy.ndarray

get_multi_index()[source]#

Returns the multi-index set of the basis.

Returns

Array of the coefficients with size (cardinality_of_basis, dimensions).

Return type

numpy.ndarray

get_parameters()[source]#

Returns the list of parameters.

Returns

Contains the Parameter objects belonging to the polynomial.

Return type

list

get_points()[source]#

Returns the samples based on the sampling strategy.

Returns

The sampled quadrature points with shape (number_of_samples, dimension).

Return type

numpy.ndarray

get_points_and_weights()[source]#

Returns the samples and weights based on the sampling strategy.

Returns

Tuple containing two ndarrays; the sampled quadrature points with shape (number_of_samples, dimension), and the corresponding quadrature weights with shape (number_of_samples, 1).

Return type

tuple

get_poly(stack_of_points, custom_multi_index=None)[source]#

Evaluates the value of each polynomial basis function at a set of points.

Parameters
  • stack_of_points (numpy.ndarray) – An ndarray with shape (number of observations, dimensions) at which the polynomial must be evaluated.

  • custom_multi_index (numpy.ndarray, optional) – Array containing a custom multi-index set, in the format given by get_elements().

Returns

Array of shape (cardinality, number_of_observations) corresponding to the polynomial basis function evaluations at the stack_of_points.

Return type

numpy.ndarray

get_poly_grad(stack_of_points, dim_index=None)[source]#

Evaluates the gradient for each of the polynomial basis functions at a set of points, with respect to each input variable.

Parameters
  • stack_of_points (numpy.ndarray) – An ndarray with shape (number_of_observations, dimensions) at which the gradient must be evaluated.

  • dim_index (int, optional) – Index of the dimension to evaluate the gradient for.

Returns

A list with d elements, where d corresponds to the dimension of the problem. Each element is a numpy.ndarray of shape (cardinality, number_of_observations) corresponding to the gradient polynomial evaluations at the stack_of_points.

Return type

list

get_poly_hess(stack_of_points)[source]#

Evaluates the Hessian for each of the polynomial basis functions at a set of points, with respect to each input variable.

Parameters

stack_of_points (numpy.ndarray) – An ndarray with shape (number_of_observations, dimensions) at which the Hessian must be evaluated.

Returns

A list with d^2 elements, where d corresponds to the dimension of the model. Each element is a numpy.ndarray of shape (cardinality, number_of_observations) corresponding to the hessian polynomial evaluations at the stack_of_points.

Return type

list

get_polyfit(stack_of_points, uq=False)[source]#

Evaluates the /polynomial approximation of a function (or model data) at prescribed points.

Parameters
  • stack_of_points (numpy.ndarray) – An ndarray with shape (number_of_observations, dimensions) at which the polynomial fit must be evaluated at.

  • uq (bool, optional) – If True, the estimated uncertainty (standard deviation) of the polynomial approximation is also returned (see [5]).

Returns

Array with shape (1, number_of_observations) corresponding to the polynomial approximation of the model.

Return type

numpy.ndarray

get_polyfit_function()[source]#

Returns a callable polynomial approximation of a function (or model data).

Returns

A callable function.

Return type

Callable

get_polyfit_grad(stack_of_points, dim_index=None)[source]#

Evaluates the gradient of the polynomial approximation of a function (or model data) at prescribed points.

Parameters
  • stack_of_points (numpy.ndarray) – An ndarray with shape (number_of_observations, dimensions) at which the polynomial fit approximation’s gradient must be evaluated at.

  • dim_index (int, optional) – Index of dimension to evaluate gradient for.

Returns

Array of shape (dimensions, number_of_observations) corresponding to the polynomial gradient approximation of the model.

Return type

numpy.ndarray

get_polyfit_grad_function()[source]#

Returns a callable for the gradients of the polynomial approximation of a function (or model data).

Returns

A callable function.

Return type

Callable

get_polyfit_hess(stack_of_points)[source]#

Evaluates the hessian of the polynomial approximation of a function (or model data) at prescribed points.

Parameters

stack_of_points (numpy.ndarray) – An ndarray with shape (number_of_observations, dimensions) at which the polynomial fit approximation’s Hessian must be evaluated at.

Returns

Array of shape (dimensions, dimensions, number_of_observations) corresponding to the polynomial Hessian approximation of the model.

Return type

numpy.ndarray

get_polyfit_hess_function()[source]#

Returns a callable for the hessian of the polynomial approximation of a function (or model data).

Returns

A callable function.

Return type

Callable

get_polyscore(X_test=None, y_test=None, metric='adjusted_r2')[source]#

Evaluates the accuracy of the polynomial approximation using the selected accuracy metric.

Training accuracy is evaluated on the data used for fitting the polynomial. Testing accuracy is evaluated on new data if it is provided by the X_test and y_test arguments (both must be provided together).

Parameters
  • X_test (numpy.ndarray, optional) – An ndarray with shape (number_of_observations, dimensions), containing new data X_test data (optional).

  • y_test (numpy.ndarray, optional) – An ndarray with shape (number_of_observations, 1) containing new y_test data (optional).

  • metric (str, optional) – An optional string containing the scoring metric to use. Avaliable options are: adjusted_r2, r2, mae, rmse, or normalised_mae (default: adjusted_r2).

Returns

  • float – If X_test and y_test are None. The training score of the model, output as a float.

  • tuple – If X_test and y_test are given. Tuple containing the train and test scores of the model.

get_skewness_and_kurtosis()[source]#

Computes the skewness and kurtosis of the model.

Returns

Tuple containing two floats; the approximated skewness and kurtosis from the polynomial fit.

Return type

tuple

get_sobol_indices(order)[source]#

Computes the Sobol’ indices.

Parameters

order (int) – The order of the Sobol’ indices required.

Returns

A dict comprising of Sobol’ indices and constitutent mixed orders of the parameters.

Return type

dict

get_summary(filename=None, tosay=False)[source]#

A simple utility that writes a file summarising what the polynomial approximation has determined.

Parameters
  • filename (str, optional) – The filename to write to. If None, the summary is written to effective-quadratures-output.txt.

  • tosay (bool, optional) – True will replace “-” signs with “minus” when writing to file for compatibility with os.say().

get_total_sobol_indices()[source]#

Computes the total Sobol’ indices.

Returns

Array of length d (number of input parameters) of total Sobol’ indices.

Return type

numpy.ndarray

get_weights()[source]#

Computes quadrature weights.

Returns

Array of the corresponding quadrature weights with shape (number_of_samples, 1).

Return type

numpy.ndarray

plot_model_vs_data(ax=None, sample_data=None, metric='adjusted_r2', show=True)[source]#

Plots the polynomial approximation against the true data. See plot_model_vs_data() for full description.

plot_parameters(ax=None, cols=2, show=True)[source]#

Plots the probability density functions for all Parameters within a Polynomial. See plot_parameters() for full description.

plot_polyfit_1D(ax=None, uncertainty=True, output_variances=None, number_of_points=200, show=True)[source]#

Plots a 1D only polynomial fit to the data. See plot_polyfit_1D() for full description.

plot_sobol(ax=None, order=1, show=True, labels=None, kwargs={})[source]#

Plots a polynomial’s Sobol’ indices of a given order. See plot_sobol() for full description.

plot_sobol_heatmap(parameters=None, show=True, ax=None)[source]#

Generates a heatmap showing the first and second order Sobol indices. See plot_sobol_heatmap() for full description.

plot_total_sobol(ax=None, show=True, labels=None, kwargs={})[source]#

Plots a polynomial’s total-order Sobol’ indices. See plot_total_sobol() for full description.

set_model(model=None, model_grads=None)[source]#

Computes the coefficients of the polynomial via the method selected.

If model evaluations and/or gradients have not yet been provided to Poly via sample-outputs and/or sample-gradients, they can be given here via the model and model_grads arguments.

Parameters
  • model (Callable,numpy.ndarray) – The function that needs to be approximated. In the absence of a callable function, the input can be a ndarray model evaluations at the quadrature points.

  • model_grads (Callable,numpy.ndarray) – The gradient of the function that needs to be approximated. In the absence of a callable gradient function, the input can be a ndarray of gradient evaluations at the quadrature points.

Examples

Defining a simple univariate example
>>> x = np.linspace(-2,2,20).reshape(-1,1)
>>> f = lambda x: x**2 - 1
>>> param = eq.Parameter(lower=-2,upper=2,order=2)
>>> basis = eq.Basis('univariate')
Least-squares regression of existing data
>>> y = f(x)
>>> poly = eq.Poly(param,basis,method='least-squares',sampling_args={'mesh': 'user-defined','sample-points':x,'sample-outputs':y})
>>> poly.set_model()
>>> eq.datasets.score(f(xtest),poly.get_polyfit(xtest),metric='rmse')
5.438959822042073e-16
Numerical integration with model provided as a callable function
>>> poly = eq.Poly(param,basis,method='numerical-integration')
>>> poly.set_model(f)
>>> eq.datasets.score(f(xtest),poly.get_polyfit(xtest),metric='rmse')
5.363652779335998e-15
Numerical integration, manually providing model evaluations at quadrature points
>>> poly = eq.Poly(param,basis,method='numerical-integration')
>>> x = poly.get_points()
>>> y = f(x)
>>> poly.set_model(y)
>>> eq.datasets.score(f(xtest),poly.get_polyfit(xtest),metric='rmse')
5.363652779335998e-15
Raises

ValueError – model values should be a column vector.

equadratures.poly.evaluate_model(points, function)[source]#

Evaluates the model function at given values.

Parameters
  • points (numpy.ndarray) – An ndarray with shape (number_of_observations, dimensions) at which the gradient must be evaluated.

  • function (Callable) – A callable argument for the function.

Returns

Array of function evaluations.

Return type

numpy.ndarray

equadratures.poly.evaluate_model_gradients(points, fungrad, format='matrix')[source]#

Evaluates the model gradient at given values.

Parameters
  • points (numpy.ndarray) – An ndarray with shape (number_of_observations, dimensions) at which the gradient must be evaluated.

  • fungrad (Callable) – A callable argument for the function’s gradients.

  • format (str, optional) – The format in which the output is to be provided: matrix will output a numpy.ndarray of shape (number_of_observations, dimensions) with gradient values, while vector will stack all the vectors in this matrix to yield a numpy.ndarray with shape (number_of_observations x dimensions, 1).

Returns

Array of gradient evaluations.

Return type

numpy.ndarray

equadratures.poly.vector_to_2D_grid(coefficients, index_set)[source]#

Handy function that converts a vector of coefficients into a matrix based on index set values.

Parameters
  • coefficients (numpy.ndarray) – An ndarray with shape (N, 1) where N corresponds to the number of coefficient values.

  • index_set (numpy.ndarray) – The multi-index set of the basis.

Returns

Tuple (x,y,z,max_order), containing the x,y values of the meshgrid, the coefficient values, and the highest order.

Return type

tuple