Nataf

Nataf is a class for transforming random variables using the Nataf transformation and calculating the correlation distortion. According to the Nataf transformation theory ([61], [62]), an n-dimensional dependent random vector \(\textbf{X}=[X_1,...,X_n]\) for which the marginal cumulative distributions \(F_i(x_i)\) and the correlattion matrix \(\textbf{C}_x=[\xi_{ij}]\) are known, can be transformed (component-wise) to standard normal random vector \(\textbf{z}=[Z_1,...,Z_n]\) with correlation matrix \(\textbf{C}_z=[\rho_{ij}]\) through the transformation:

\[Z_{i}=\Phi^{-1}(F_i(X_{i}))\]

where \(\Phi(\cdot)\) is the standard normal cumulative distribution function.

This transformation causes a correlation distortion; the correlation coefficient between the standard normal variables \(Z_i\) and \(Z_j\), denoted \(\rho_{ij}\), is not equal to its counterpart in the parameter space (\(\rho_{ij} \neq \xi_{ij}\)).

If the Gaussian correlation \(\rho_{ij}\) is know, the non-Gaussian correlation \(\xi_{ij}\) can be determined through the following integral equation:

\[\xi_{ij} =\dfrac{1}{\sigma_{X_i}\sigma_{X_j}}\int_{-\infty}^{\infty}\int_{-\infty}^{\infty}\left(X_i-\mu_{X_i}\right)\left(X_j-\mu_{X_j}\right)\phi_2(Z_i,Z_j;\rho_{ij})dZ_idZ_j\]

where \(X_i =F_i^{-1}(\Phi(Z_{i}))\) and \(\phi_2(\cdot)\) is the bivariate standard normal probability density function. The integration is directly evaluated using a quadratic two-dimensional Gauss-Legendre integration. However, in the case where the non-Gaussian correlation is known \(\xi_{ij}\), the integral above cannot be inverted to solve for the Gaussian correlation \(\rho_{ij}\). In such case, iterative methods must be employed such as the Iterative Translation Approximation Method (ITAM) [55], used in UQpy.

The Jacobian of the transformation can be also estimated with the UQpy class as:

\[\mathbf{J_{xz}} = \dfrac{\partial\mathbf{x}}{\partial\mathbf{z}} =\left[\dfrac{\phi(Z_i)}{f_i(X_i)}\right]\mathbf{H}.\]

where \(\textbf{H}\) is the lower diagonal matrix resulting from the Cholesky decomposition of the correlation matrix (\(\mathbf{C_Z}\)).

The Nataf class also allows for the inverse of the Nataf transformation, i.e. transforming a vector of standard normal vector \(\textbf{z}=[Z_1,...,Z_n]\) to random variables with prescribed marginal cumulative distributions and correlation matrix \(\textbf{C}_x=[\rho_{ij}]\) according to:

\[X_{i}=F_i^{-1}(\Phi(Z_{i}))\]

Nataf Class

The Nataf class is imported using the following command:

>>> from UQpy.transformations.Nataf import Nataf

Methods

class Nataf(distributions, samples_x=None, samples_z=None, jacobian=False, corr_z=None, corr_x=None, itam_beta=1.0, itam_threshold1=0.001, itam_threshold2=0.1, itam_max_iter=100, n_gauss_points=128)[source]

Transform random variables using the Nataf or Inverse Nataf transformation

Parameters:
  • distributions (Union[Distribution, list[Distribution]]) – Probability distribution of each random variable. Must be an object of type DistributionContinuous1D or JointIndependent.

  • samples_x (Optional[ndarray]) – Random vector of shape (nsamples, n_dimensions) with prescribed probability distributions. If samples_x is provided, the Nataf class transforms them to samples_z.

  • samples_z (Optional[ndarray]) – Standard normal random vector of shape (nsamples_, n_dimensions) If samples_z is provided, the Nataf class transforms them to samples_x.

  • jacobian (bool) – A boolean whether to return the jacobian of the transformation. Default: False

  • corr_z (Optional[ndarray]) – The correlation matrix (\(\mathbf{C_Z}\)) of the standard normal random vector Z . Default: The identity matrix. If corr_z is specified, the Nataf class computes the correlation distortion integral above to solve for corr_x.

  • corr_x (Optional[ndarray]) – The correlation matrix (\(\mathbf{C_X}\)) of the random vector X . Default: The identity matrix. If corr_x is specified, the Nataf class invokes the ITAM to compute corr_z.

  • itam_beta (Union[float, int]) – A parameter selected to optimize convergence speed and desired accuracy of the ITAM method. Default: \(1.0\)

  • itam_threshold1 (Union[float, int]) – If corr_x is specified, this specifies the threshold value for the error in the ITAM method (see utilities module) given by \(\epsilon_1 = ||\mathbf{C_X}^{target} - \mathbf{C_X}^{computed}||\) Default: \(0.001\)

  • itam_threshold2 (Union[float, int]) – If corr_x is specified, this specifies the threshold value for the error difference between iterations in the ITAM method (see utilities module) given by \(\epsilon_1^{i} - \epsilon_1^{i-1}\) for iteration \(i\). Default: \(0.01\)

  • itam_max_iter (int) – Maximum number of iterations for the ITAM method. Default: \(100\)

  • n_gauss_points (int) – The number of integration points used for the numerical integration of the correlation matrix (\(\mathbf{C_Z}\)) of the standard normal random vector Z

run(samples_x=None, samples_z=None, jacobian=False)[source]

Execute the Nataf transformation or its inverse. If samples_x is provided, the run() method performs the Nataf transformation. If samples_z is provided, the run() method performs the inverse Nataf transformation.

Parameters:
  • samples_x (Optional[ndarray]) – Random vector X with prescribed probability distributions or standard normal random vector Z of shape (nsamples, n_dimensions).

  • samples_z (Optional[ndarray]) – Standard normal random vector of shape (nsamples, n_dimensions) If samples_z is provided, the Nataf class transforms them to samples_x.

  • jacobian (bool) – The jacobian of the transformation of shape (n_dimensions, n_dimensions). Default: False

static itam(distributions, corr_x, itam_max_iter=100, itam_beta=1.0, itam_threshold1=0.001, itam_threshold2=0.01)[source]

Calculate the correlation matrix \(\mathbf{C_Z}\) of the standard normal random vector \(\mathbf{Z}\) given the correlation matrix \(\mathbf{C_X}\) of the random vector \(\mathbf{X}\) using the ITAM method [55].

Parameters:
  • distributions (Union[DistributionContinuous1D, JointIndependent, list[Union[DistributionContinuous1D, JointIndependent]]]) – Probability distribution of each random variable. Must be an object of type DistributionContinuous1D or JointIndependent. distributions must have a cdf method.

  • corr_x – The correlation matrix (\(\mathbf{C_X}\)) of the random vector X . Default: The identity matrix.

  • itam_max_iter (int) – Maximum number of iterations for the ITAM method. Default: \(100\)

  • itam_beta (Union[float, int]) – A parameters selected to optimize convergence speed and desired accuracy of the ITAM method (see [62]). Default: \(1.0\)

  • itam_threshold1 (Union[float, int]) – A threshold value for the relative difference between the non-Gaussian correlation function and the underlying Gaussian. Default: \(0.001\)

  • itam_threshold2 (Union[float, int]) – If corr_x is specified, this specifies the threshold value for the error difference between iterations in the ITAM method (see utilities module) given by: \(\epsilon_1^{i} - \epsilon_1^{i-1}\) for iteration \(i\). Default: \(0.01\)

Returns:

Distorted correlation matrix (\(\mathbf{C_Z}\)) of the standard normal vector Z, List of ITAM errors for each iteration, List of ITAM difference errors for each iteration

static distortion_z2x(distributions, corr_z, n_gauss_points=1024)[source]

This is a method to calculate the correlation matrix \(\mathbf{C_x}\) of the random vector \(\mathbf{x}\) given the correlation matrix \(\mathbf{C_z}\) of the standard normal random vector \(\mathbf{z}\).

Parameters:
  • distributions (Union[Distribution, list[Distribution]]) – This is a method to calculate the correlation matrix \(\mathbf{C_x}\) of the random vector \(\mathbf{x}\) given the correlation matrix \(\mathbf{C_z}\) of the standard normal random vector \(\mathbf{z}\). This method is part of the Nataf class.

  • corr_z (ndarray) – The correlation matrix (\(\mathbf{C_z}\)) of the standard normal vector Z . Default: The identity matrix.

  • n_gauss_points (int) – The number of integration points used for the numerical integration of the correlation matrix (\(\mathbf{C_Z}\)) of the standard normal random vector Z

Returns:

Distorted correlation matrix (\(\mathbf{C_X}\)) of the random vector X.

rvs(nsamples)[source]

Generate realizations from the joint pdf of the random vector X.

Parameters:

nsamples (int) – Number of samples to generate.

Returns:

Random vector in the parameter space of shape (nsamples, n_dimensions).

Attributes

Nataf.samples_x: ndarray

Random vector of shape (nsamples, n_dimensions) with prescribed probability distributions.

Nataf.samples_z: ndarray

Standard normal random vector of shape (nsamples, n_dimensions)

Nataf.jzx: ndarray

The Jacobian of the transformation of shape (n_dimensions, n_dimensions).

Nataf.jxz: ndarray

The Jacobian of the transformation of shape (n_dimensions, n_dimensions).

Nataf.corr_z: ndarray

Distorted correlation matrix (\(\mathbf{C_Z}\)) of the standard normal vector Z.

Nataf.corr_x: ndarray

Distorted correlation matrix (\(\mathbf{C_X}\)) of the random vector X.

Nataf.H: ndarray

The lower triangular matrix resulting from the Cholesky decomposition of the correlation matrix \(\mathbf{C_Z}\).

Examples