Transformations¶
This module contains functionality for isoprobabilistic transformations in UQpy.
The module currently contains the following classes:
Nataf: Class to perform the Nataf isoprobabilistic transformations.Correlate: Class to induce correlation to a standard normal vector.Decorrelate: Class to remove correlation from a standard normal vector.
Nataf¶
Nataf is a class for transforming random variables using the Nataf transformation and calculating the correlation distortion. According to the Nataf transformation theory (1, 2), a 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:
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:
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) 3, used in UQpy.
The Jacobian of the transformation can be also estimated with the Nataf class as:
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:
Nataf Class Descriptions¶
-
class
UQpy.Transformations.Nataf(dist_object, samples_x=None, samples_z=None, jacobian=False, corr_z=None, corr_x=None, beta=None, itam_threshold1=None, itam_threshold2=None, itam_max_iter=None, verbose=False)[source]¶ Transform random variables using the Nataf or Inverse Nataf transformation
Inputs:
- dist_object ((list of )
Distributionobject(s)): Probability distribution of each random variable. Must be an object of type
DistributionContinuous1DorJointInd.
- dist_object ((list of )
- corr_x (ndarray):
The correlation matrix (\(\mathbf{C_X}\)) of the random vector X .
Default: The identity matrix.
If
corr_xis specified, theNatafclass invokes the ITAM to computecorr_z.
- corr_z (ndarray):
The correlation matrix (\(\mathbf{C_Z}\)) of the standard normal random vector Z .
Default: The identity matrix.
If
corr_zis specified, theNatafclass computes the correlation distortion integral above to solve forcorr_x.
- itam_threshold1 (float):
If
corr_xis specified, this specifies the threshold value for the error in the ITAM method (seeUtilitiesmodule) given by:\[\epsilon_1 = ||\mathbf{C_X}^{target} - \mathbf{C_X}^{computed}||\]Default: 0.001
- itam_threshold2 (float):
If
corr_xis specified, this specifies the threshold value for the error difference between iterations in the ITAM method (seeUtilitiesmodule) given by:\[\epsilon_1^{i} - \epsilon_1^{i-1}\]for iteration \(i\).
Default: 0.01
- beta (float):
A parameter selected to optimize convergence speed and desired accuracy of the ITAM method.
Default: 1.0
- samples_x (ndarray):
Random vector of shape
(nsamples, dimension)with prescribed probability distributions.If samples_x is provided, the
Natafclass transforms them to samples_z.
- samples_z (ndarray):
Standard normal random vector of shape
(nsamples, dimension)If samples_z is provided, the
Natafclass transforms them to samples_x.
- jacobian (‘Boolean’):
A boolean whether to return the jacobian of the transformation.
Default: False
- itam_max_iter (int):
Maximum number of iterations for the ITAM method.
Default: 100
- verbose (Boolean):
A boolean declaring whether to write text to the terminal.
Default:
False
Attributes:
- corr_z (ndarray):
Distorted correlation matrix (\(\mathbf{C_Z}\)) of the standard normal vector Z.
- corr_x (ndarray):
Distorted correlation matrix (\(\mathbf{C_X}\)) of the random vector X.
- H (ndarray):
The lower triangular matrix resulting from the Cholesky decomposition of the correlation matrix \(\mathbf{C_Z}\).
- itam_error1 (list)
List of ITAM errors for each iteration
- itam_error2 (list)
List of ITAM difference errors for each iteration
- samples_x (ndarray):
Random vector of shape
(nsamples, dimension)with prescribed probability distributions.
- samples_z (ndarray):
Standard normal random vector of shape
(nsamples, dimension)
- jxz (ndarray):
The Jacobian of the transformation of shape
(dimension, dimension).
- jzx (ndarray):
The Jacobian of the transformation of shape
(dimension, dimension).
Methods:
-
static
distortion_z2x(dist_object, corr_z, verbose=None)[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}\).
This method is part of the
Natafclass.Inputs:
- dist_object ((list of )
Distributionobject(s)): Probability distribution of each random variable. Must be an object of type
DistributionContinuous1DorJointInd.dist_object must have an
icdfmethod.
- dist_object ((list of )
- corr_z (ndarray):
The correlation matrix (\(\mathbf{C_z}\)) of the standard normal vector Z .
Default: The
identitymatrix.
- verbose (Boolean):
A boolean declaring whether to write text to the terminal.
Default:
False
Output/Returns:
- corr_x (ndarray):
Distorted correlation matrix (\(\mathbf{C_X}\)) of the random vector X.
-
static
itam(dist_object, corr_x, itam_max_iter=None, beta=None, itam_threshold1=None, itam_threshold2=None, verbose=None)[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 3.
Inputs:
- dist_object ((list of )
Distributionobject(s)): Probability distribution of each random variable. Must be an object of type
DistributionContinuous1DorJointInd.dist_object must have a
cdfmethod.
- dist_object ((list of )
- corr_x (ndarray):
The correlation matrix (\(\mathbf{C_X}\)) of the random vector X .
Default: The
identitymatrix.
- itam_threshold1 (float):
If
corr_xis specified, this specifies the threshold value for the error in the ITAM method (seeUtilitiesmodule) given by:\[\epsilon_1 = ||\mathbf{C_X}^{target} - \mathbf{C_X}^{computed}||\]
- itam_max_iter (int):
Maximum number of iterations for the ITAM method.
Default: 100
- itam_threshold1 (float):
A threshold value for the relative difference between the non-Gaussian correlation function and the underlying Gaussian.
Default: 0.001
- itam_threshold2 (float):
If
corr_xis specified, this specifies the threshold value for the error difference between iterations in the ITAM method (seeUtilitiesmodule) given by:\[\epsilon_1^{i} - \epsilon_1^{i-1}\]for iteration \(i\).
Default: 0.01
- beta (float):
A parameters selected to optimize convergence speed and desired accuracy of the ITAM method (see 2).
Default: 1.0
- verbose (Boolean):
A boolean declaring whether to write text to the terminal.
Default: False
Output/Returns:
- corr_z (ndarray):
Distorted correlation matrix (\(\mathbf{C_Z}\)) of the standard normal vector Z.
- itam_error1 (list)
List of ITAM errors for each iteration
- itam_error2 (list)
List of ITAM difference errors for each iteration
-
run(samples_x=None, samples_z=None, jacobian=False)[source]¶ Execute the Nataf transformation or its inverse.
If samples_x is provided, the
runmethod performs the Nataf transformation. If samples_z is provided, therunmethod performs the inverse Nataf transformation.** Input:**
- samples_x or samples_z (ndarray):
Random vector X with prescribed probability distributions or standard normal random vector Z of shape``(nsamples, dimension)``.
- jacobian (Boolean):
The jacobian of the transformation of shape
(dimension, dimension).Default:
False
Correlate¶
Correlate is a class to induce correlation to an uncorrelated standard normal vector \(\textbf{u}=[U_1,...,U_n]\), given the correlation matrix \(\textbf{C}_z=[\rho_{ij}]\). The correlated standard normal vector \(\textbf{z}=[Z_1,...,Z_n]\) can be calculated as:
where \(\mathbf{H}\) is the lower triangular matrix resulting from the Cholesky decomposition of the correlation matrix, i.e. \(\mathbf{C_z}=\mathbf{H}\mathbf{H}^\intercal\).
Correlate Class Descriptions¶
-
class
UQpy.Transformations.Correlate(samples_u, corr_z)[source]¶ A class to induce correlation to standard normal random variables.
Inputs:
- samples_u (ndarray):
Uncorrelated standard normal vector of shape
(nsamples, dimension).
- corr_z (ndarray):
The correlation matrix (\(\mathbf{C_Z}\)) of the standard normal random vector Z .
Attributes:
- samples_z (ndarray):
Correlated standard normal vector of shape
(nsamples, dimension).
- H (ndarray):
The lower diagonal matrix resulting from the Cholesky decomposition of the correlation matrix (\(\mathbf{C_Z}\)).
Decorrelate¶
Decorrelate is a class to remove correlation from an correlated standard normal vector \(\textbf{z}=[Z_1,...,Z_n]\) with correlation matrix \(\textbf{C}_z=[\rho_{ij}]\). The uncorrelated standard normal vector \(\textbf{u}=[U_1,...,U_n]\) can be calculated as:
Decorrelate Class Descriptions¶
-
class
UQpy.Transformations.Decorrelate(samples_z, corr_z)[source]¶ A class to remove correlation from correlated standard normal random variables.
Inputs:
- samples_z (ndarray):
Correlated standard normal vector of shape
(nsamples, dimension).
- corr_z (ndarray):
The correlation matrix (\(\mathbf{C_Z}\)) of the standard normal random vector Z .
Attributes:
- samples_u (ndarray):
Uncorrelated standard normal vector of shape
(nsamples, dimension).
- H (ndarray):
The lower diagonal matrix resulting from the Cholesky decomposition of the correlation matrix (\(\mathbf{C_Z}\)).
References:
- 1
A. Nataf, Determination des distributions dont les marges sont donnees, C. R. Acad. Sci. vol. 225, pp. 42-43, Paris, 1962.
- 2(1,2)
R. Lebrun and A. Dutfoy, An innovating analysis of the Nataf transformation from the copula viewpoint, Prob. Eng. Mech., vol. 24, pp. 312-320, 2009.
- 3(1,2)
Kim, H. and Shields, M.D. (2015). “Modeling strongly non-Gaussian non-stationary stochastic processes using the Iterative Translation Approximation Method and Karhunen-Loeve Expansion,” Computers and Structures. 161: 31-42.