Matrix

class Matrix(*args)

Real rectangular matrix.

Parameters:
n_rowsint, n_r > 0, optional

Number of rows. Default is 1.

n_columnsint, n_c > 0, optional

Number of columns. Default is 1.

valuessequence of float with size n_r \times n_c, optional

Values. column-major ordering is used (like Fortran) for reshaping the flat list of values. Default creates a zero matrix.

Examples

Create a matrix

>>> import openturns as ot
>>> M = ot.Matrix(2, 2, range(2 * 2))
>>> print(M)
[[ 0 2 ]
 [ 1 3 ]]

Get or set terms

>>> print(M[0, 0])
0.0
>>> M[0, 0] = 1.
>>> print(M[0, 0])
1.0
>>> print(M[:, 0])
[[ 1 ]
 [ 1 ]]

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

>>> import numpy as np
>>> np_2d_array = np.array([[1.0, 2.0, 3.0], [4.0, 5.0, 6.0]])
>>> ot_matrix = ot.Matrix(np_2d_array)

and back

>>> np_matrix = np.matrix(ot_matrix)

Basic linear algebra operations (provided the dimensions are compatible)

>>> A = ot.Matrix([[1.0, 2.0], [3.0, 4.0], [5.0, 6.0]])
>>> B = ot.Matrix(np.eye(2))
>>> C = ot.Matrix(3, 2, [1.] * 3 * 2)
>>> print(A * B - C)
[[ 0 1 ]
 [ 2 3 ]
 [ 4 5 ]]

Methods

clean(threshold)

Set elements smaller than a threshold to zero.

computeGram([transpose])

Compute the associated Gram matrix.

computeHadamardProduct(other)

Compute the Hadamard product matrix.

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.

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.

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.

isEmpty()

Tell if the matrix is empty.

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

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

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 ]]
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
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]
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
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
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 rectangular linear system whose the present matrix is the operator.

Parameters:
rhsPoint 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 rectangular linear system.

Notes

This will handle both matrices and vectors, as well as underdetermined, square or overdetermined linear systems although you’d better type explicitly your matrix if it has some properties that could simplify the resolution (see TriangularMatrix, SquareMatrix).

This uses LAPACK’s DGELSY. The RCOND parameter of this routine can be changed through the MatrixImplementation-DefaultSmallPivot key of the ResourceMap.

Examples

>>> import openturns as ot
>>> import numpy as np
>>> M = ot.Matrix([[1.0, 2.0], [3.0, 4.0], [5.0, 6.0]])
>>> b = ot.Point([1.0] * 3)
>>> 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:
MTMatrix

The transposed matrix.

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.transpose())
[[ 1 3 5 ]
 [ 2 4 6 ]]

Examples using the class

Estimate Wilks and empirical quantile

Estimate Wilks and empirical quantile

Test independence

Test independence

Fit a parametric copula

Fit a parametric copula

Fit a non parametric copula

Fit a non parametric copula

Estimate a multivariate ARMA process

Estimate a multivariate ARMA process

Estimate a non stationary covariance function

Estimate a non stationary covariance function

Visualize pairs

Visualize pairs

Visualize clouds

Visualize clouds

Create and draw multivariate distributions

Create and draw multivariate distributions

Quick start guide to distributions

Quick start guide to distributions

Draw minimum volume level sets

Draw minimum volume level sets

Create a copula

Create a copula

Extract the copula from a distribution

Extract the copula from a distribution

Assemble copulas

Assemble copulas

Create a functional basis process

Create a functional basis process

Create a parametric spectral density function

Create a parametric spectral density function

Export a field to VTK

Export a field to VTK

Create a stationary covariance model

Create a stationary covariance model

Create a custom covariance model

Create a custom covariance model

Create a discrete Markov chain process

Create a discrete Markov chain process

Create a stationary covariance model

Create a stationary covariance model

Create a normal process

Create a normal process

Draw a field

Draw a field

Sample trajectories from a Gaussian Process with correlated outputs

Sample trajectories from a Gaussian Process with correlated outputs

Draw fields

Draw fields

Mixture of experts

Mixture of experts

Kriging :configure the optimization solver

Kriging :configure the optimization solver

Metamodel of a field function

Metamodel of a field function

Estimate moments from Taylor expansions

Estimate moments from Taylor expansions

Use the post-analytical importance sampling algorithm

Use the post-analytical importance sampling algorithm

Subset Sampling

Subset Sampling

Time variant system reliability problem

Time variant system reliability problem

Create unions or intersections of events

Create unions or intersections of events

An illustrated example of a FORM probability estimate

An illustrated example of a FORM probability estimate

Create an event based on a process

Create an event based on a process

Estimate Sobol indices on a field to point function

Estimate Sobol indices on a field to point function

Parallel coordinates graph as sensitivity tool

Parallel coordinates graph as sensitivity tool

Use the ANCOVA indices

Use the ANCOVA indices

Create mixed deterministic and probabilistic designs of experiments

Create mixed deterministic and probabilistic designs of experiments

Create a quadratic function

Create a quadratic function

Increase the input dimension of a function

Increase the input dimension of a function

Function manipulation

Function manipulation

Calibration without observed inputs

Calibration without observed inputs

Calibration of the deflection of a tube

Calibration of the deflection of a tube

Calibration of the flooding model

Calibration of the flooding model

Calibration of the Chaboche mechanical model

Calibration of the Chaboche mechanical model

Bayesian calibration of a computer code

Bayesian calibration of a computer code

Bayesian calibration of the flooding model

Bayesian calibration of the flooding model

Linear Regression with interval-censored observations

Linear Regression with interval-censored observations

Iterated Functions System

Iterated Functions System

Compute leave-one-out error of a polynomial chaos expansion

Compute leave-one-out error of a polynomial chaos expansion

Compute confidence intervals of a regression model from data

Compute confidence intervals of a regression model from data

Compute confidence intervals of a univariate noisy function

Compute confidence intervals of a univariate noisy function

Control algorithm termination

Control algorithm termination

A quick start guide to graphs

A quick start guide to graphs