CovarianceMatrix

class CovarianceMatrix(*args)

Covariance (real symmetric positive definite) matrix.

Parameters
sizeint, n > 0, optional

Matrix size. Default is 1.

valuessequence of float with size n^2, optional

Values. OpenTURNS uses column-major ordering (like Fortran) for reshaping the flat list of values. Default creates an identity matrix.

Raises
TypeErrorIf the matrix is not symmetric.

Examples

Create a matrix

>>> import openturns as ot
>>> C = ot.CovarianceMatrix(2, [1.0, 0.5, 0.5, 1.0])
>>> print(C)
[[ 1   0.5 ]
 [ 0.5 1   ]]

Get or set terms

>>> print(C[0, 1])
0.5
>>> C[0, 1] = 0.6
>>> print(C[0, 1])
0.6
>>> print(C[:, 0])
[[ 1   ]
 [ 0.6 ]]

Create an openturns matrix from a symmetric numpy 2d-array (or matrix, or 2d-list)…

>>> import numpy as np
>>> np_2d_array = np.array([[1.0, 0.5], [0.5, 1.0]])
>>> ot_matrix = ot.CovarianceMatrix(np_2d_array)

and back

>>> np_matrix = np.matrix(ot_matrix)

Methods

checkSymmetry()

Check if the internal representation is really symmetric.

clean(threshold)

Set elements smaller than a threshold to zero.

computeCholesky([keepIntact])

Compute the Cholesky factor.

computeDeterminant([keepIntact])

Compute the determinant.

computeEV([keepIntact])

Compute the eigenvalues decomposition (EVD).

computeEigenValues([keepIntact])

Compute eigenvalues.

computeGram([transpose])

Compute the associated Gram matrix.

computeLogAbsoluteDeterminant([keepIntact])

Compute the logarithm of the absolute value of the determinant.

computeQR([fullQR, keepIntact])

Compute the QR factorization.

computeSVD([fullSVD, keepIntact])

Compute the singular values decomposition (SVD).

computeSingularValues([keepIntact])

Compute the singular values.

computeTrace()

Compute the trace of the matrix.

getClassName()

Accessor to the object’s name.

getDimension()

Accessor to the dimension (the number of rows).

getId()

Accessor to the object’s id.

getImplementation()

Accessor to the underlying implementation.

getName()

Accessor to the object’s name.

getNbColumns()

Accessor to the number of columns.

getNbRows()

Accessor to the number of rows.

isDiagonal()

Test whether the matrix is diagonal or not.

isEmpty()

Tell if the matrix is empty.

isPositiveDefinite()

Test whether the matrix is positive definite or not.

reshape(newRowDim, newColDim)

Reshape the matrix.

reshapeInPlace(newRowDim, newColDim)

Reshape the matrix, in place.

setName(name)

Accessor to the object’s name.

solveLinearSystem(*args)

Solve a square linear system whose the present matrix is the operator.

transpose()

Transpose the matrix.

__init__(*args)

Initialize self. See help(type(self)) for accurate signature.

checkSymmetry()

Check if the internal representation is really symmetric.

clean(threshold)

Set elements smaller than a threshold to zero.

Parameters
thresholdfloat

Threshold for zeroing elements.

Returns
cleaned_matrixMatrix

Input matrix with elements smaller than the threshold set to zero.

computeCholesky(keepIntact=True)

Compute the Cholesky factor.

The Cholesky factor of a covariance (real symmetric positive definite) matrix \mat{C} is the lower triangular matrix \mat{L} such that:

\mat{C} = \mat{L} \Tr{\mat{L}}

Parameters
keep_intactbool, optional

A flag telling whether the present matrix can be overwritten or not. Default is True and leaves the present matrix unchanged.

Returns
cholesky_factorSquareMatrix

The left (lower) Cholesky factor.

Notes

This uses LAPACK’s DPOTRF.

computeDeterminant(keepIntact=True)

Compute the determinant.

Parameters
keep_intactbool, optional

A flag telling whether the present matrix can be overwritten or not. Default is True and leaves the present matrix unchanged.

Returns
determinantfloat

The square matrix determinant.

Examples

>>> import openturns as ot
>>> A = ot.SquareMatrix([[1.0, 2.0], [3.0, 4.0]])
>>> A.computeDeterminant()
-2.0
computeEV(keepIntact=True)

Compute the eigenvalues decomposition (EVD).

The eigenvalues decomposition of a square matrix \mat{M} with size n reads:

\mat{M} = \mat{\Phi} \mat{\Lambda} \Tr{\mat{\Phi}}

where \mat{\Lambda} is an n \times n diagonal matrix and \mat{\Phi} is an n \times n orthogonal matrix.

Parameters
keep_intactbool, optional

A flag telling whether the present matrix can be overwritten or not. Default is True and leaves the present matrix unchanged.

Returns
eigenvaluesPoint

The vector of eigenvalues with size n that form the diagonal of the n \times n matrix \mat{\Lambda} of the EVD.

PhiSquareComplexMatrix

The left matrix of the EVD.

Notes

This uses LAPACK’S DSYEV.

Examples

>>> import openturns as ot
>>> import numpy as np
>>> M = ot.SymmetricMatrix([[1.0, 2.0], [2.0, -4.0]])
>>> eigen_values, Phi = M.computeEV()
>>> Lambda = ot.SquareMatrix(M.getDimension())
>>> for i in range(eigen_values.getSize()):
...     Lambda[i, i] = eigen_values[i]
>>> np.testing.assert_array_almost_equal(Phi * Lambda * Phi.transpose(), M)
computeEigenValues(keepIntact=True)

Compute eigenvalues.

Parameters
keep_intactbool, optional

A flag telling whether the present matrix can be overwritten or not. Default is True and leaves the present matrix unchanged.

Returns
eigenvaluesPoint

Eigenvalues.

See also

computeEV

Examples

>>> import openturns as ot
>>> M = ot.SymmetricMatrix([[1.0, 2.0], [2.0, -4.0]])
>>> print(M.computeEigenValues())
[-4.70156,1.70156]
computeGram(transpose=True)

Compute the associated Gram matrix.

Parameters
transposedbool

Tells if matrix is to be transposed or not. Default value is True

Returns
MMTMatrix

The Gram matrix.

Notes

When transposed is set to True, the method computes cM^t \times \cM. Otherwise it computes \cM \ times \cM^t

Examples

>>> import openturns as ot
>>> M = ot.Matrix([[1.0, 2.0], [3.0, 4.0], [5.0, 6.0]])
>>> MtM = M.computeGram()
>>> print(MtM)
[[ 35 44 ]
 [ 44 56 ]]
>>> MMt = M.computeGram(False)
>>> print(MMt)
[[  5 11 17 ]
 [ 11 25 39 ]
 [ 17 39 61 ]]
computeLogAbsoluteDeterminant(keepIntact=True)

Compute the logarithm of the absolute value of the determinant.

Parameters
keep_intactbool, optional

A flag telling whether the present matrix can be overwritten or not. Default is True and leaves the present matrix unchanged.

Returns
determinantfloat

The logarithm of the absolute value of the square matrix determinant.

signfloat

The sign of the determinant.

Examples

>>> import openturns as ot
>>> A = ot.SquareMatrix([[1.0, 2.0], [3.0, 4.0]])
>>> A.computeLogAbsoluteDeterminant()
[0.693147..., -1.0]
computeQR(fullQR=False, keepIntact=True)

Compute the QR factorization. By default, it is the economic decomposition which is computed.

The economic QR factorization of a rectangular matrix \mat{M} with n_r \geq n_c (more rows than columns) is defined as follows:

\mat{M} = \mat{Q} \mat{R}
        = \mat{Q} \begin{bmatrix} \mat{R_1} \\ \mat{0} \end{bmatrix}
        = \begin{bmatrix} \mat{Q_1}, \mat{Q_2} \end{bmatrix}
          \begin{bmatrix} \mat{R_1} \\ \mat{0} \end{bmatrix}
        = \mat{Q_1} \mat{R_1}

where \mat{R_1} is an n_c \times n_c upper triangular matrix, \mat{Q_1} is n_r \times n_c, \mat{Q_2} is n_r \times (n_r - n_c), and \mat{Q_1} and \mat{Q_2} both have orthogonal columns.

Parameters
full_qrbool, optional

A flag telling whether Q, R or Q1, R1 are returned. Default is False and returns Q1, R1.

keep_intactbool, optional

A flag telling whether the present matrix is preserved or not in the computation of the decomposition. Default is True and leaves the present matrix unchanged.

Returns
Q1Matrix

The orthogonal matrix of the economic QR factorization.

R1TriangularMatrix

The right (upper) triangular matrix of the economic QR factorization.

QMatrix

The orthogonal matrix of the full QR factorization.

RTriangularMatrix

The right (upper) triangular matrix of the full QR factorization.

Notes

The economic QR factorization is often used for solving overdetermined linear systems (where the operator \mat{M} has n_r \geq n_c) in the least-square sense because it implies solving a (simple) triangular system:

\vect{\hat{x}} = \arg\min\limits_{\vect{x} \in \Rset^{n_r}} \|\mat{M} \vect{x} - \vect{b}\|
               = \mat{R_1}^{-1} (\Tr{\mat{Q_1}} \vect{b})

This uses LAPACK’s DGEQRF and DORGQR.

Examples

>>> import openturns as ot
>>> import numpy as np
>>> M = ot.Matrix([[1.0, 2.0], [3.0, 4.0], [5.0, 6.0]])
>>> Q1, R1 = M.computeQR()
>>> np.testing.assert_array_almost_equal(Q1 * R1, M)
computeSVD(fullSVD=False, keepIntact=True)

Compute the singular values decomposition (SVD).

The singular values decomposition of a rectangular matrix \mat{M} with size n_r > n_c reads:

\mat{M} = \mat{U} \mat{\Sigma} \Tr{\mat{V}}

where \mat{U} is an n_r \times n_r orthogonal matrix, \mat{\Sigma} is an n_r \times n_c diagonal matrix and \mat{V} is an n_c \times n_c orthogonal matrix.

Parameters
fullSVDbool, optional

Whether the null parts of the orthogonal factors are explicitely stored or not. Default is False and computes a reduced SVD.

keep_intactbool, optional

A flag telling whether the present matrix can be overwritten or not. Default is True and leaves the present matrix unchanged.

Returns
singular_valuesPoint

The vector of singular values with size n = \min(n_r, n_c) that form the diagonal of the n_r \times n_c matrix \mat{\Sigma} of the SVD.

USquareMatrix

The left orthogonal matrix of the SVD.

VTSquareMatrix

The transposed right orthogonal matrix of the SVD.

Notes

This uses LAPACK’s DGESDD.

Examples

>>> import openturns as ot
>>> import numpy as np
>>> M = ot.Matrix([[1.0, 2.0], [3.0, 4.0], [5.0, 6.0]])
>>> singular_values, U, VT = M.computeSVD(True)
>>> Sigma = ot.Matrix(M.getNbRows(), M.getNbColumns())
>>> for i in range(singular_values.getSize()):
...     Sigma[i, i] = singular_values[i]
>>> np.testing.assert_array_almost_equal(U * Sigma * VT, M)
computeSingularValues(keepIntact=True)

Compute the singular values.

Parameters
fullSVDbool, optional

Whether the null parts of the orthogonal factors are explicitely stored or not. Default is False and computes a reduced SVD.

keep_intactbool, optional

A flag telling whether the present matrix can be overwritten or not. Default is True and leaves the present matrix unchanged.

Returns
singular_valuesPoint

The vector of singular values with size n = \min(n_r, n_c) that form the diagonal of the n_r \times n_c matrix \mat{\Sigma} of the SVD decomposition.

See also

computeSVD

Examples

>>> import openturns as ot
>>> M = ot.Matrix([[1.0, 2.0], [3.0, 4.0], [5.0, 6.0]])
>>> print(M.computeSingularValues(True))
[9.52552,0.514301]
computeTrace()

Compute the trace of the matrix.

Returns
tracefloat

The trace of the matrix.

Examples

>>> import openturns as ot
>>> M = ot.SquareMatrix([[1.0, 2.0], [3.0, 4.0]])
>>> M.computeTrace()
5.0
getClassName()

Accessor to the object’s name.

Returns
class_namestr

The object class name (object.__class__.__name__).

getDimension()

Accessor to the dimension (the number of rows).

Returns
dimensionint
getId()

Accessor to the object’s id.

Returns
idint

Internal unique identifier.

getImplementation()

Accessor to the underlying implementation.

Returns
implImplementation

The implementation class.

getName()

Accessor to the object’s name.

Returns
namestr

The name of the object.

getNbColumns()

Accessor to the number of columns.

Returns
n_columnsint
getNbRows()

Accessor to the number of rows.

Returns
n_rowsint
isDiagonal()

Test whether the matrix is diagonal or not.

Returns
testbool

Answer.

isEmpty()

Tell if the matrix is empty.

Returns
is_emptybool

True if the matrix contains no element.

Examples

>>> import openturns as ot
>>> M = ot.Matrix([[]])
>>> M.isEmpty()
True
isPositiveDefinite()

Test whether the matrix is positive definite or not.

A matrix \mat{M} is positive definite if \Tr{\vect{z}} \mat{M} \vect{z} is positive for every compatible non-zero column vector \vect{z}.

Returns
testbool

Answer.

Notes

This uses LAPACK’s DPOTRF.

reshape(newRowDim, newColDim)

Reshape the matrix.

Parameters
newRowDimint

The row dimension of the reshaped matrix.

newColDimint

The column dimension of the reshaped matrix.

Returns
MTMatrix

The reshaped matrix.

Notes

If the size of the reshaped matrix is smaller than the size of the matrix to be reshaped, only the newRowDim\times newColDim first elements are kept (in a column-major storage sense). If the size is greater, the new elements are set to zero.

Examples

>>> import openturns as ot
>>> M = ot.Matrix([[1.0, 2.0], [3.0, 4.0], [5.0, 6.0]])
>>> print(M)
[[ 1 2 ]
 [ 3 4 ]
 [ 5 6 ]]
>>> print(M.reshape(1, 6))
1x6
[[ 1 3 5 2 4 6 ]]
>>> print(M.reshape(2, 2))
[[ 1 5 ]
 [ 3 2 ]]
>>> print(M.reshape(2, 6))
2x6
[[ 1 5 4 0 0 0 ]
 [ 3 2 6 0 0 0 ]]
reshapeInPlace(newRowDim, newColDim)

Reshape the matrix, in place.

Parameters
newRowDimint

The row dimension of the reshaped matrix.

newColDimint

The column dimension of the reshaped matrix.

Notes

If the size of the reshaped matrix is smaller than the size of the matrix to be reshaped, only the newRowDim\times newColDim first elements are kept (in a column-major storage sense). If the size is greater, the new elements are set to zero. If the size is unchanged, no copy of data is done.

Examples

>>> import openturns as ot
>>> M = ot.Matrix([[1.0, 2.0], [3.0, 4.0], [5.0, 6.0]])
>>> print(M)
[[ 1 2 ]
 [ 3 4 ]
 [ 5 6 ]]
>>> M.reshapeInPlace(1, 6)
>>> print(M)
1x6
[[ 1 3 5 2 4 6 ]]
>>> M.reshapeInPlace(2, 2)
>>> print(M)
[[ 1 5 ]
 [ 3 2 ]]
>>> M.reshapeInPlace(2, 6)
>>> print(M)
2x6
[[ 1 5 0 0 0 0 ]
 [ 3 2 0 0 0 0 ]]
setName(name)

Accessor to the object’s name.

Parameters
namestr

The name of the object.

solveLinearSystem(*args)

Solve a square linear system whose the present matrix is the operator.

Parameters
rhssequence of float or Matrix with n_r values or rows, respectively

The right hand side member of the linear system.

keep_intactbool, optional

A flag telling whether the present matrix can be overwritten or not. Default is True and leaves the present matrix unchanged.

Returns
solutionPoint or Matrix

The solution of the square linear system.

Notes

This will handle both matrices and vectors. Note that you’d better type explicitely the matrix if it has some properties that could simplify the resolution (see TriangularMatrix).

This uses LAPACK’S DGESV for matrices and DGELSY for vectors.

Examples

>>> import openturns as ot
>>> import numpy as np
>>> M = ot.SquareMatrix([[1.0, 2.0], [3.0, 4.0]])
>>> b = ot.Point([1.0] * 2)
>>> x = M.solveLinearSystem(b)
>>> np.testing.assert_array_almost_equal(M * x, b)
transpose()

Transpose the matrix.

Returns
MTSquareMatrix

The transposed matrix.

Examples

>>> import openturns as ot
>>> M = ot.SquareMatrix([[1.0, 2.0], [3.0, 4.0]])
>>> print(M)
[[ 1 2 ]
 [ 3 4 ]]
>>> print(M.transpose())
[[ 1 3 ]
 [ 2 4 ]]