Correlations#

Utilities for dealing with correlated inputs.

class equadratures.correlations.Correlations(correlation_matrix, poly=None, parameters=None, method=None, verbose=False)[source]#

The class defines methods for polynomial approximations with correlated inputs, including the Nataf transform and Gram-Schmidt process.

Parameters
  • correlation_matrix (numpy.ndarray) – The correlation matrix associated with the joint distribution.

  • poly (Poly, optional) – Polynomial defined with parameters with marginal distributions in uncorrelated space.

  • parameters (list, optional) – List of parameters with marginal distributions.

  • method (str, optional) – nataf-transform or gram-schmidt.

  • verbose (bool, optional) – Display Cholesky decomposition of the fictive matrix.

Example

>>> def func(x):
>>>     s1 = s[0]
>>>     s2 = s[1]
>>>     return s1**2 - s2**3 - 2. * s1 - 0.5 * s2
>>>
>>> parameters = [Parameter(distribution='beta', shape_parameter_A=5.0, shape_parameter_B=2.0, lower=0.0, upper=1.0, order=15) for _ in range(2)]
>>> basis = Basis('total-order')
>>> uncorr_poly = Poly(parameters, basis, method='least-squares',
>>>                     sampling_args={'mesh': 'monte-carlo',
>>>                                     'subsampling-algorithm': 'lu'})
>>>
>>> corr_mat = np.array([[1.0, -0.855],
>>>                     [-0.855, 1.0]])
>>>
>>> corr_gs = Correlations(corr_mat, poly=uncorr_poly, method='gram-schmidt')
>>> corr_gs.set_model(func)
>>> corrected_poly_gs = corr_gs.get_transformed_poly()
>>> print(corrected_poly_gs.get_mean_and_variance())

References

  1. Melchers, R. E., (1945) Structural Reliability Analysis and Predictions. John Wiley and Sons, second edition.

  2. Jakeman, J. D. et al., (2019) Polynomial chaos expansions for dependent random variables.

get_correlated_samples(N=None, X=None)[source]#

Method for generating correlated samples.

Parameters
  • N (int) – Number of correlated samples required.

  • X (numpy.ndarray, optional) – Points in the uncorrelated space to map to the correlated space.

Returns

Array with shape (N, M), which contains the correlated samples.

Return type

numpy.ndarray

get_pdf(X)[source]#

Evaluate PDF at the sample points.

Parameters

X (numpy.ndarray) – Array of sample points with shape (number_of_points, dimensions).

Returns

Array with shape (N,) containing evaluations of the PDF.

Return type

numpy.ndarray

get_points()[source]#

Returns the quadrature points accounting for correlations.

For Nataf transform, returns points in the correlated standard normal space. For Gram-Schmidt, returns points according to the GS polynomial basis.

Returns

Array of quadrature points with shape (number_of_samples, dimension).

Return type

numpy.ndarray

get_transformed_poly()[source]#

Returns the transformed polynomial.

Returns

The transformed polynomial.

Return type

Poly

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

Computes the coefficients of transformed polynomial (equivalent to calling self.get_transformed_poly().set_model(...))

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

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