TriangularMatrix

class TriangularMatrix(*args)

Hermitian Matrix.

Available constructors:

TriangularMatrix(dim)

TriangularMatrix(dim, isLower)

Parameters:
diminteger

The dimension of the triangular matrix (square matrix with dim rows and dim columns).

isLowerbool

Flag telling if the matrix is triangular lower (True) or upper (False). Default is True.

See also

Matrix

Notes

The triangular matrix is filled with 0.

Methods

clean(threshold)

Set elements smaller than a threshold to zero.

computeDeterminant([keepIntact])

Compute the determinant.

computeEV([keepIntact])

Compute the eigenvalues decomposition (EVD).

computeEigenValues([keepIntact])

Compute eigenvalues.

computeGram([transpose])

Compute the associated Gram matrix.

computeHadamardProduct(other)

Compute the Hadamard product matrix.

computeLargestEigenValueModule(*args)

Compute the largest eigenvalue module.

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.

computeSumElements()

Compute the sum of the matrix elements.

computeTrace()

Compute the trace of the matrix.

getClassName()

Accessor to the object's name.

getDiagonal([k])

Get the k-th diagonal of the matrix.

getDiagonalAsPoint([k])

Get the k-th diagonal of the matrix.

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.

isLowerTriangular()

Test whether the matrix is lower triangular or upper triangular.

reshape(newRowDim, newColDim)

Reshape the matrix.

reshapeInPlace(newRowDim, newColDim)

Reshape the matrix, in place.

setDiagonal(*args)

Set the k-th diagonal of the matrix.

setName(name)

Accessor to the object's name.

solveLinearSystem(*args)

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

solveLinearSystemInPlace(*args)

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

squareElements()

Square the Matrix, ie each element of the matrix is squared.

transpose()

Transpose the matrix.

checkTriangularity

__init__(*args)
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.

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} \mat{\Phi}^{-1}

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:
eigen_valuesComplexCollection

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 DGEEV.

Examples

>>> import openturns as ot
>>> import numpy as np
>>> M = ot.SquareMatrix([[1.0, 2.0], [3.0, 4.0]])
>>> eigen_values, Phi = M.computeEV()
>>> Lambda = ot.SquareComplexMatrix(M.getDimension())
>>> for i in range(eigen_values.getSize()):
...     Lambda[i, i] = eigen_values[i]
>>> # from scipy.linalg import inv # SquareComplexMatrix does not implement solveLinearSystem
>>> # Phi, Lambda = np.matrix(Phi), np.matrix(Lambda)
>>> # np.testing.assert_array_almost_equal(Phi * Lambda * inv(Phi), 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:
eigenvaluesComplexCollection

Eigenvalues.

See also

computeEV

Examples

>>> import openturns as ot
>>> M = ot.SquareMatrix([[1.0, 2.0], [3.0, 4.0]])
>>> M.computeEigenValues()
[(-0.372281,0),(5.37228,0)]
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 ]]
computeHadamardProduct(other)

Compute the Hadamard product matrix.

Parameters:
MatMatrix

The right hand matrix

Returns:
CMatrix

The Hadamard product matrix.

Notes

The matrix \cC resulting from the Hadamard product of the matrices \cA and :cB is as follows:

\cC_{i,j} = \cA_{i,j} * \cB_{i,j}

Examples

>>> import openturns as ot
>>> A = ot.Matrix([[1.0, 2.0], [3.0, 4.0]])
>>> B = ot.Matrix([[1.0, 2.0], [3.0, 4.0]])
>>> C = A.computeHadamardProduct(B)
>>> print(C)
[[  1  4 ]
 [  9 16 ]]
>>> print(B.computeHadamardProduct(A))
[[  1  4 ]
 [  9 16 ]]
computeLargestEigenValueModule(*args)

Compute the largest eigenvalue module.

Parameters:
maximumIterationsint, optional

The maximum number of power iterations to perform to get the approximation. Default is given by the ‘Matrix-LargestEigenValueIterations’ key in the ResourceMap.

epsilonfloat, optional

The target relative error. Default is given by the ‘Matrix-LargestEigenValueRelativeError’ key in the ResourceMap.

Returns:
largestEigenvalueModulefloat

The largest eigenvalue module.

Examples

>>> import openturns as ot
>>> M = ot.SquareMatrix([[1.0, 2.0], [3.0, 4.0]])
>>> M.computeLargestEigenValueModule()
5.3722...
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 explicitly 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 explicitly 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]
computeSumElements()

Compute the sum of the matrix elements.

Returns:
suma float

The sum of the elements.

Notes

We compute here the sum of elements of the matrix \cM, that defines as:

s = \sum_{i=1}^{nRows}\sum_{j=1}^{nColumns} \cM_{i,j}

Examples

>>> import openturns as ot
>>> M = ot.Matrix([[1.0, 2.0], [3.0, 4.0], [5.0, 6.0]])
>>> s = M.computeSumElements()
>>> print(s)
21.0
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__).

getDiagonal(k=0)

Get the k-th diagonal of the matrix.

Parameters:
kint

The k-th diagonal to extract Default value is 0

Returns:
D: Matrix

The k-th diagonal.

Examples

>>> import openturns as ot
>>> M = ot.Matrix([[1.0, 2.0, 3.0], [4.0, 5.0, 6.0], [7.0, 8.0, 9.0]])
>>> diag = M.getDiagonal()
>>> print(diag)
[[ 1 ]
 [ 5 ]
 [ 9 ]]
>>> print(M.getDiagonal(1))
[[ 2 ]
 [ 6 ]]
getDiagonalAsPoint(k=0)

Get the k-th diagonal of the matrix.

Parameters:
kint

The k-th diagonal to extract Default value is 0

Returns:
ptPoint

The k-th digonal.

Examples

>>> import openturns as ot
>>> M = ot.Matrix([[1.0, 2.0, 3.0], [4.0, 5.0, 6.0], [7.0, 8.0, 9.0]])
>>> pt = M.getDiagonalAsPoint()
>>> print(pt)
[1,5,9]
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

A copy of the underlying implementation object.

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
isLowerTriangular()

Test whether the matrix is lower triangular or upper triangular.

Returns:
isLowerbool

Flag telling if the matrix is triangular lower (True) or upper (False).

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 ]]
setDiagonal(*args)

Set the k-th diagonal of the matrix.

Parameters:
diagPoint or Matrix or a float

Value(s) used to fill the diagonal of the matrix

kint

The k-th diagonal to fill Default value is 0

Examples

>>> import openturns as ot
>>> M = ot.Matrix([[1.0, 2.0, 3.0], [4.0, 5.0, 6.0], [7.0, 8.0, 9.0]])
>>> M.setDiagonal([-1, 23, 9])
>>> print(M)
[[ -1  2  3 ]
 [  4 23  6 ]
 [  7  8  9 ]]
>>> M.setDiagonal(1.0)
>>> print(M)
[[ 1 2 3 ]
 [ 4 1 6 ]
 [ 7 8 1 ]]
>>> M.setDiagonal([2, 6, 9])
>>> print(M)
[[ 2 2 3 ]
 [ 4 6 6 ]
 [ 7 8 9 ]]
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.

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 explicitly 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 = [1.0] * 2
>>> x = M.solveLinearSystem(b)
>>> np.testing.assert_array_almost_equal(M * x, b)
solveLinearSystemInPlace(*args)

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

Similar to solveLinearSystem() except the matrix is modified in-place during the resolution avoiding the need to allocate an extra copy if the original copy is not re-used.

squareElements()

Square the Matrix, ie each element of the matrix is squared.

Examples

>>> import openturns as ot
>>> M = ot.Matrix([[1.0, 2.0], [3.0, 4.0], [5.0, 6.0]])
>>> M.squareElements()
>>> print(M)
[[  1  4 ]
 [  9 16 ]
 [ 25 36 ]]
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 ]]