Matrix¶
- class Matrix(*args)¶
Real rectangular matrix.
- Parameters:
- n_rowsint, , optional
Number of rows. Default is 1.
- n_columnsint, , optional
Number of columns. Default is 1.
- valuessequence of float with size , 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.
Compute the sum of the matrix elements.
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.
Accessor to the underlying implementation.
getName
()Accessor to the object's name.
Accessor to the number of columns.
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.
Square the Matrix, ie each element of the matrix is squared.
Transpose the matrix.
- __init__(*args)¶
- clean(threshold)¶
Set elements smaller than a threshold to zero.
- Parameters:
- thresholdfloat
Threshold for zeroing elements.
- Returns:
- cleaned_matrix
Matrix
Input matrix with elements smaller than the threshold set to zero.
- cleaned_matrix
- computeGram(transpose=True)¶
Compute the associated Gram matrix.
- Parameters:
- transposedbool
Tells if matrix is to be transposed or not. Default value is True
- Returns:
- MMT
Matrix
The Gram matrix.
- MMT
Notes
When transposed is set to True, the method computes . Otherwise it computes
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.
Notes
The matrix resulting from the Hadamard product of the matrices and :cB is as follows:
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 with (more rows than columns) is defined as follows:
where is an upper triangular matrix, is , is , and and 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:
- Q1
Matrix
The orthogonal matrix of the economic QR factorization.
- R1
TriangularMatrix
The right (upper) triangular matrix of the economic QR factorization.
- Q
Matrix
The orthogonal matrix of the full QR factorization.
- R
TriangularMatrix
The right (upper) triangular matrix of the full QR factorization.
- Q1
Notes
The economic QR factorization is often used for solving overdetermined linear systems (where the operator has ) in the least-square sense because it implies solving a (simple) triangular system:
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 with size reads:
where is an orthogonal matrix, is an diagonal matrix and is an 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_values
Point
The vector of singular values with size that form the diagonal of the matrix of the SVD.
- U
SquareMatrix
The left orthogonal matrix of the SVD.
- VT
SquareMatrix
The transposed right orthogonal matrix of the SVD.
- singular_values
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_values
Point
The vector of singular values with size that form the diagonal of the matrix of the SVD decomposition.
- singular_values
See also
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 , that defines as:
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.
- D:
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:
- pt
Point
The k-th digonal.
- pt
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:
- MT
Matrix
The reshaped matrix.
- MT
Notes
If the size of the reshaped matrix is smaller than the size of the matrix to be reshaped, only the 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 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:
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:
- Returns:
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 ]]
Examples using the class¶
Estimate Wilks and empirical quantile
Estimate a multivariate ARMA process
Estimate a non stationary covariance function
Create and draw multivariate distributions
Quick start guide to distributions
Draw minimum volume level sets
Extract the copula from a distribution
Create a functional basis process
Create a parametric spectral density function
Create a stationary covariance model
Create a custom covariance model
Create a discrete Markov chain process
Create a stationary covariance model
Sample trajectories from a Gaussian Process with correlated outputs
Kriging :configure the optimization solver
Estimate moments from Taylor expansions
Use the post-analytical importance sampling algorithm
Time variant system reliability problem
Create unions or intersections of events
An illustrated example of a FORM probability estimate
Create an event based on a process
Estimate Sobol indices on a field to point function
Parallel coordinates graph as sensitivity tool
Create mixed deterministic and probabilistic designs of experiments
Increase the input dimension of a function
Calibration without observed inputs
Calibration of the deflection of a tube
Calibration of the flooding model
Calibration of the Chaboche mechanical model
Bayesian calibration of a computer code
Bayesian calibration of the flooding model
Linear Regression with interval-censored observations
Compute confidence intervals of a regression model from data
Compute confidence intervals of a univariate noisy function