Parameter#

Definition of a univariate parameter.

class equadratures.parameter.Parameter(order=1, distribution='Uniform', endpoints=None, shape_parameter_A=None, shape_parameter_B=None, variable='parameter', lower=None, upper=None, weight_function=None)[source]#

This class defines a univariate parameter.

Parameters
  • lower (float, optional) – Lower bound for the parameter.

  • upper (float, optional) – Upper bound for the parameter.

  • order (int, optional) – Order of the parameter.

  • param_type (str, optional) – The type of distribution that characterizes the parameter (see [1, 2]). Options include chebyshev (arcsine), gaussian, truncated-gaussian, beta, cauchy, exponential, uniform, triangular, gamma, weibull, rayleigh, pareto, lognormal, students-t, logistic, gumbel, chi and chi-squared. If no string is provided, a uniform distribution is assumed. Data-driven and custom analytical parameters can also be constructed by setting this option to data and analytical and providing a weight_function (see examples).

  • shape_parameter_A (float, optional) – Most of the aforementioned distributions are characterized by two shape parameters. For instance, in the case of a gaussian (or truncated-gaussian), this represents the mean. In the case of a beta distribution this represents the alpha value. For a uniform distribution this input is not required.

  • shape_parameter_B (float, optional) – This is the second shape parameter that characterizes the distribution selected. In the case of a gaussian or truncated-gaussian, this is the variance.

  • data (numpy.ndarray, optional) – A data-set with shape (number_of_data_points, 2), where the first column comprises of parameter values, while the second column corresponds to the data observations. This input should only be used with the Analytical distribution.

  • endpoints (str, optional) – If set to both, then the quadrature points and weights will have end-points, based on Gauss-Lobatto quadrature rules. If set to upper or lower a Gauss-Radau rule is used to compute one end-point at either the upper or lower bound.

  • weight_function (Weight, optional) – An instance of Weight, which contains a bespoke analytical or data-driven weight (probability density) function.

Examples

A uniform parameter
>>> param = eq.Parameter(distribution='uniform', lower=-2, upper=2., order=3)
A beta parameter
>>> param = eq.Parameter(distribution='beta', lower=-2., upper=15., order=4,
>>>        shape_parameter_A=3.2, shape_parameter_B=1.7)
A data-driven parameter
>>> pdf = eq.Weight( stats.gaussian_kde(data, bw_method='silverman'),
>>>        support=[-3, 3.2])
>>> param = eq.Parameter(distribution='analytical',
>>>        weight_function=pdf, order=2)

References

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

  2. Gautschi, W., (1985) Orthogonal Polynomials-Constructive Theory and Applications. Journal of Computational and Applied Mathematics 12 (1985), pp. 61-76. Paper

get_cdf(points=None)[source]#

Computes the cumulative density function associated with the Parameter.

Parameters

points (numpy.ndarray, optional) – Values of the parameter at which the CDF must be evaluated.

Returns

  • numpy.ndarray – If points!=None. ndarray containing the cumulative density function evaluated at the points in points.

  • tuple – If points=None. A tuple (x, cdf), where cdf is the cumulative density function evaluated at the points in x.

get_description()[source]#

Provides a description of the Parameter.

Returns

A description of the parameter.

Return type

str

get_icdf(cdf_values)[source]#

Computes the inverse cumulative density function associated with the Parameter.

Parameters

cdf_values (numpy.ndarray) – Values of the cumulative density function for which its inverse needs to be computed.

Returns

The inverse cumulative density function.

Return type

numpy.ndarray

get_jacobi_eigenvectors(order=None)[source]#

Computes the eigenvectors of the Jacobi matrix.

Parameters

order (int) – Order of the recurrence coefficients.

Returns

Array of eigenvectors.

Return type

numpy.ndarray

get_jacobi_matrix(order=None, ab=None)[source]#

Computes the Jacobi matrix—a tridiagonal matrix of the recurrence coefficients.

Parameters

order (int) – Order of the recurrence coefficients.

Returns

2D array containing the Jacobi matrix.

Return type

numpy.ndarray

get_pdf(points=None)[source]#

Computes the probability density function associated with the Parameter.

Parameters

points (numpy.ndarray, optional) – Values of the parameter at which the PDF must be evaluated.

Returns

  • numpy.ndarray – If points!=None. ndarray containing the probability density function evaluated at the points in points.

  • tuple – If points=None. A tuple (x, pdf), where pdf is the probability density function evaluated at the points in x.

get_recurrence_coefficients(order=None)[source]#

Generates the recurrence coefficients.

Parameters

order (int, optional) – Order of the recurrence coefficients.

Returns

Array of recurrence coefficients.

Return type

numpy.ndarray

get_samples(number_of_samples_required)[source]#

Generates samples from the distribution associated with the Parameter.

Parameters

number_of_samples_required (int) – Number of samples that are required.

Returns

The generated samples.

Return type

numpy.ndarray

plot_cdf(ax=None, show=True, lim_range=True)[source]#

Plots the cumulative density function for a Parameter. See plot_cdf() for full description.

plot_orthogonal_polynomials(ax=None, order_limit=None, number_of_points=200, show=True)[source]#

Plots the first few orthogonal polynomials. See plot_orthogonal_polynomials() for full description.

plot_pdf(ax=None, data=None, show=True, lim_range=True)[source]#

Plots the probability density function for a Parameter. See plot_pdf() for full description.