Arbitrary precision complex ball matrices using Arb#
AUTHORS:
Clemens Heuberger (2014-10-25): Initial version.
This is a rudimentary binding to the Arb library; it may be useful to refer to its documentation for more details.
- class sage.matrix.matrix_complex_ball_dense.Matrix_complex_ball_dense#
Bases:
sage.matrix.matrix_dense.Matrix_dense
Matrix over a complex ball field. Implemented using the
acb_mat
type of the Arb library.EXAMPLES:
sage: MatrixSpace(CBF, 3)(2) [2.000000000000000 0 0] [ 0 2.000000000000000 0] [ 0 0 2.000000000000000] sage: matrix(CBF, 1, 3, [1, 2, -3]) [ 1.000000000000000 2.000000000000000 -3.000000000000000]
- charpoly(var='x', algorithm=None)#
Compute the characteristic polynomial of this matrix.
EXAMPLES:
sage: from sage.matrix.benchmark import hilbert_matrix sage: mat = hilbert_matrix(5).change_ring(ComplexBallField(10)) sage: mat.charpoly() x^5 + ([-1.8 +/- 0.0258])*x^4 + ([0.3 +/- 0.05...)*x^3 + ([+/- 0.0...])*x^2 + ([+/- 0.0...])*x + [+/- 0.0...]
- contains(other)#
Test if the set of complex matrices represented by
self
is contained in that represented byother
.EXAMPLES:
sage: b = CBF(0, RBF(0, rad=.1r)); b [+/- 0.101]*I sage: matrix(CBF, [0, b]).contains(matrix(CBF, [0, 0])) True sage: matrix(CBF, [0, b]).contains(matrix(CBF, [b, 0])) False sage: matrix(CBF, [b, b]).contains(matrix(CBF, [b, 0])) True
- determinant()#
Compute the determinant of this matrix.
EXAMPLES:
sage: matrix(CBF, [[1/2, 1/3], [1, 1]]).determinant() [0.1666666666666667 +/- ...e-17] sage: matrix(CBF, [[1/2, 1/3], [1, 1]]).det() [0.1666666666666667 +/- ...e-17] sage: matrix(CBF, [[1/2, 1/3]]).determinant() Traceback (most recent call last): ... ValueError: self must be a square matrix
- eigenvalues(other=None, extend=None)#
(Experimental.) Compute rigorous enclosures of the eigenvalues of this matrix.
INPUT:
self
– an \(n \times n\) matrixother
– unsupported (generalized eigenvalue problem), should beNone
extend
– ignored
OUTPUT:
A
Sequence
of complex balls of length equal to the size of the matrix.Each element represents one eigenvalue with the correct multiplicities in case of overlap. The output intervals are either disjoint or identical, and identical intervals are guaranteed to be grouped consecutively. Each complete run of \(k\) identical balls thus represents a cluster of exactly \(k\) eigenvalues which could not be separated from each other at the current precision, but which could be isolated from the other eigenvalues.
There is currently no guarantee that the algorithm converges as the working precision is increased.
See the Arb documentation for more information.
EXAMPLES:
sage: from sage.matrix.benchmark import hilbert_matrix sage: mat = hilbert_matrix(5).change_ring(CBF) sage: mat.eigenvalues() doctest:...: FutureWarning: This class/method/function is marked as experimental. ... [[1.567050691098...] + [+/- ...]*I, [0.208534218611...] + [+/- ...]*I, [3.287928...e-6...] + [+/- ...]*I, [0.000305898040...] + [+/- ...]*I, [0.011407491623...] + [+/- ...]*I] sage: mat = Permutation([2, 1, 4, 5, 3]).to_matrix().dense_matrix().change_ring(CBF) sage: mat.eigenvalues() Traceback (most recent call last): ... ValueError: unable to certify the eigenvalues sage: precond = matrix(ZZ, [[-1, -2, 2, 2, -2], [2, -2, -2, -2, 2], ....: [-2, 2, -1, 2, 1], [2, 1, -1, 0, 2], [-2, 0, 1, -1, 1]]) sage: (~precond*mat*precond).eigenvalues() [[-0.5000000000000...] + [-0.8660254037844...]*I, [-1.000000000000...] + [+/- ...]*I, [-0.5000000000000...] + [0.8660254037844...]*I, [1.000000000000...] + [+/- ...]*I, [1.000000000000...] + [+/- ...]*I]
See also
- eigenvectors_left(other=None, extend=True)#
(Experimental.) Compute rigorous enclosures of the eigenvalues and left eigenvectors of this matrix.
INPUT:
self
– an \(n \times n\) matrixother
– unsupported (generalized eigenvalue problem), should beNone
extend
– ignored
OUTPUT:
A list of triples of the form
(eigenvalue, [eigenvector], 1)
.Unlike
eigenvalues()
andeigenvectors_left_approx()
, this method currently fails in the presence of multiple eigenvalues.Additionally, there is currently no guarantee that the algorithm converges as the working precision is increased.
See the Arb documentation for more information.
EXAMPLES:
sage: mat = matrix(CBF, 3, [2, 3, 5, 7, 11, 13, 17, 19, 23]) sage: eigval, eigvec, _ = mat.eigenvectors_left()[0] sage: eigval [1.1052996349...] + [+/- ...]*I sage: eigvec[0] ([0.69817246751...] + [+/- ...]*I, [-0.67419514369...] + [+/- ...]*I, [0.240865343781...] + [+/- ...]*I) sage: eigvec[0] * (mat - eigval) ([+/- ...] + [+/- ...]*I, [+/- ...] + [+/- ...]*I, [+/- ...] + [+/- ...]*I)
- eigenvectors_left_approx(other=None, extend=None)#
(Experimental.) Compute non-rigorous approximations of the left eigenvalues and eigenvectors of this matrix.
INPUT:
self
– an \(n \times n\) matrixother
– unsupported (generalized eigenvalue problem), should beNone
extend
– ignored
OUTPUT:
A list of triples of the form
(eigenvalue, [eigenvector], 1)
. The eigenvalue and the entries of the eigenvector are complex balls with zero radius.No guarantees are made about the accuracy of the output.
See the Arb documentation for more information.
EXAMPLES:
sage: mat = matrix(CBF, 3, [2, 3, 5, 7, 11, 13, 17, 19, 23]) sage: eigval, eigvec, _ = mat.eigenvectors_left_approx()[0] sage: eigval [1.1052996349... +/- ...] sage: eigvec[0] ([0.69817246751...], [-0.67419514369...], [0.240865343781...]) sage: eigvec[0] * (mat - eigval) ([+/- ...], [+/- ...], [+/- ...])
See also
- eigenvectors_right(other=None, extend=None)#
(Experimental.) Compute rigorous enclosures of the eigenvalues and eigenvectors of this matrix.
INPUT:
self
– an \(n \times n\) matrixother
– unsupported (generalized eigenvalue problem), should beNone
extend
– ignored
OUTPUT:
A list of triples of the form
(eigenvalue, [eigenvector], 1)
.Unlike
eigenvalues()
andeigenvectors_right_approx()
, this method currently fails in the presence of multiple eigenvalues.Additionally, there is currently no guarantee that the algorithm converges as the working precision is increased.
See the Arb documentation for more information.
EXAMPLES:
sage: from sage.matrix.benchmark import hilbert_matrix sage: mat = hilbert_matrix(3).change_ring(CBF) sage: eigval, eigvec, _ = mat.eigenvectors_right()[0] doctest:...: FutureWarning: This class/method/function is marked as experimental. ... sage: eigval [1.40831892712...] + [+/- ...]*I sage: eigvec [([0.82704492697...] + [+/- ...]*I, [0.45986390436...] + [+/- ...]*I, [0.32329843524...] + [+/- ...]*I)] sage: (mat - eigval)*eigvec[0] ([+/- ...] + [+/- ...]*I, [+/- ...] + [+/- ...]*I, [+/- ...] + [+/- ...]*I)
See also
- eigenvectors_right_approx(other=None, extend=None)#
(Experimental.) Compute non-rigorous approximations of the eigenvalues and eigenvectors of this matrix.
INPUT:
self
– an \(n \times n\) matrixother
– unsupported (generalized eigenvalue problem), should beNone
extend
– ignored
OUTPUT:
A list of triples of the form
(eigenvalue, [eigenvector], 1)
. The eigenvalue and the entries of the eigenvector are complex balls with zero radius.No guarantees are made about the accuracy of the output.
See the Arb documentation for more information.
EXAMPLES:
sage: from sage.matrix.benchmark import hilbert_matrix sage: mat = hilbert_matrix(3).change_ring(CBF) sage: eigval, eigvec, _ = mat.eigenvectors_right_approx()[0] doctest:...: FutureWarning: This class/method/function is marked as experimental. ... sage: eigval [1.40831892712...] sage: eigval.rad() 0.00000000 sage: eigvec [([0.8270449269720...], [0.4598639043655...], [0.3232984352444...])] sage: (mat - eigval)*eigvec[0] ([1e-15 +/- ...], [2e-15 +/- ...], [+/- ...])
See also
- exp()#
Compute the exponential of this matrix.
EXAMPLES:
sage: matrix(CBF, [[i*pi, 1], [0, i*pi]]).exp() [[-1.00000000000000 +/- ...e-16] + [+/- ...e-16]*I [-1.00000000000000 +/- ...e-16] + [+/- ...e-16]*I] [ 0 [-1.00000000000000 +/- ...e-16] + [+/- ...e-16]*I] sage: matrix(CBF, [[1/2, 1/3]]).exp() Traceback (most recent call last): ... ValueError: self must be a square matrix
- identical(other)#
Test if the corresponding entries of two complex ball matrices represent the same balls.
EXAMPLES:
sage: a = matrix(CBF, [[1/3,2],[3,4]]) sage: b = matrix(CBF, [[1/3,2],[3,4]]) sage: a == b False sage: a.identical(b) True
- overlaps(other)#
Test if two matrices with complex ball entries represent overlapping sets of complex matrices.
EXAMPLES:
sage: b = CBF(0, RBF(0, rad=0.1r)); b [+/- 0.101]*I sage: matrix(CBF, [0, b]).overlaps(matrix(CBF, [b, 0])) True sage: matrix(CBF, [1, 0]).overlaps(matrix(CBF, [b, 0])) False
- trace()#
Compute the trace of this matrix.
EXAMPLES:
sage: matrix(CBF, [[1/3, 1/3], [1, 1]]).trace() [1.333333333333333 +/- ...e-16] sage: matrix(CBF, [[1/2, 1/3]]).trace() Traceback (most recent call last): ... ValueError: self must be a square matrix
- transpose()#
Return the transpose of
self
.EXAMPLES:
sage: m = matrix(CBF, 2, 3, [1, 2, 3, 4, 5, 6]) sage: m.transpose() [1.000000000000000 4.000000000000000] [2.000000000000000 5.000000000000000] [3.000000000000000 6.000000000000000] sage: m.transpose().parent() Full MatrixSpace of 3 by 2 dense matrices over Complex ball field with 53 bits of precision