Dense matrices over the Complex Double Field using NumPy¶
EXAMPLES:
sage: b = Mat(CDF,2,3).basis()
sage: b[0,0]
[1.0 0.0 0.0]
[0.0 0.0 0.0]
>>> from sage.all import *
>>> b = Mat(CDF,Integer(2),Integer(3)).basis()
>>> b[Integer(0),Integer(0)]
[1.0 0.0 0.0]
[0.0 0.0 0.0]
We deal with the case of zero rows or zero columns:
sage: m = MatrixSpace(CDF,0,3)
sage: m.zero_matrix()
[]
>>> from sage.all import *
>>> m = MatrixSpace(CDF,Integer(0),Integer(3))
>>> m.zero_matrix()
[]
AUTHORS:
Jason Grout (2008-09): switch to NumPy backend
Josh Kantor
William Stein: many bug fixes and touch ups.
- class sage.matrix.matrix_complex_double_dense.Matrix_complex_double_dense[source]¶
Bases:
Matrix_double_dense
Class that implements matrices over the real double field. These are supposed to be fast matrix operations using C doubles. Most operations are implemented using numpy which will call the underlying BLAS on the system.
EXAMPLES:
sage: # needs sage.symbolic sage: m = Matrix(CDF, [[1,2*I],[3+I,4]]) sage: m**2 [-1.0 + 6.0*I 10.0*I] [15.0 + 5.0*I 14.0 + 6.0*I] sage: n= m^(-1); n # abs tol 1e-15 [ 0.3333333333333333 + 0.3333333333333333*I 0.16666666666666669 - 0.16666666666666666*I] [-0.16666666666666666 - 0.3333333333333333*I 0.08333333333333331 + 0.08333333333333333*I]
>>> from sage.all import * >>> # needs sage.symbolic >>> m = Matrix(CDF, [[Integer(1),Integer(2)*I],[Integer(3)+I,Integer(4)]]) >>> m**Integer(2) [-1.0 + 6.0*I 10.0*I] [15.0 + 5.0*I 14.0 + 6.0*I] >>> n= m**(-Integer(1)); n # abs tol 1e-15 [ 0.3333333333333333 + 0.3333333333333333*I 0.16666666666666669 - 0.16666666666666666*I] [-0.16666666666666666 - 0.3333333333333333*I 0.08333333333333331 + 0.08333333333333333*I]
To compute eigenvalues, use the methods
left_eigenvectors()
orright_eigenvectors()
:sage: p,e = m.right_eigenvectors() # needs sage.symbolic
>>> from sage.all import * >>> p,e = m.right_eigenvectors() # needs sage.symbolic
The result is a pair
(p,e)
, wherep
is a diagonal matrix of eigenvalues ande
is a matrix whose columns are the eigenvectors.To solve a linear system \(Ax = b\) where
A = [[1,2*I],[3+I,4]]
andb = [5,6]
:sage: b = vector(CDF,[5,6]) # needs sage.symbolic sage: m.solve_right(b) # abs tol 1e-14 # needs sage.symbolic (2.6666666666666665 + 0.6666666666666669*I, -0.3333333333333333 - 1.1666666666666667*I)
>>> from sage.all import * >>> b = vector(CDF,[Integer(5),Integer(6)]) # needs sage.symbolic >>> m.solve_right(b) # abs tol 1e-14 # needs sage.symbolic (2.6666666666666665 + 0.6666666666666669*I, -0.3333333333333333 - 1.1666666666666667*I)
See the methods
QR()
,LU()
, andSVD()
for QR, LU, and singular value decomposition.