Schur {Matrix} | R Documentation |
Computes the Schur decomposition and eigenvalues of a square matrix; see the BACKGROUND information below.
Schur(x, vectors, ...)
x |
numeric
square Matrix (inheriting from class |
vectors |
logical. When |
... |
further arguments passed to or from other methods. |
Based on the Lapack subroutine dgees
.
If vectors
are TRUE
, as per default:
If x
is a Matrix
an object of class Schur
, otherwise, for a
traditional matrix
x
, a list
with
components T
, Q
, and EValues
.
If vectors
are FALSE
, a list with components
T |
the upper quasi-triangular (square) matrix of the Schur decomposition. |
EValues |
If A
is a square matrix, then A = Q T t(Q)
, where
Q
is orthogonal, and T
is upper block-triangular
(nearly triangular with either 1 by 1 or 2 by 2 blocks on the
diagonal) where the 2 by 2 blocks correspond to (non-real) complex
eigenvalues.
The eigenvalues of A
are the same as those of T
,
which are easy to compute. The Schur form is used most often for
computing non-symmetric eigenvalue decompositions, and for computing
functions of matrices such as matrix exponentials.
Anderson, E., et al. (1994). LAPACK User's Guide, 2nd edition, SIAM, Philadelphia.
Schur(Hilbert(9)) # Schur factorization (real eigenvalues) (A <- Matrix(round(rnorm(5*5, sd = 100)), nrow = 5)) (Sch.A <- Schur(A)) eTA <- eigen(Sch.A@T) str(SchA <- Schur(A, vectors=FALSE))# no 'T' ==> simple list stopifnot(all.equal(eTA$values, eigen(A)$values, tolerance = 1e-13), all.equal(eTA$values, local({z <- Sch.A@EValues z[order(Mod(z), decreasing=TRUE)]}), tolerance = 1e-13), identical(SchA$T, Sch.A@T), identical(SchA$EValues, Sch.A@EValues)) ## For the faint of heart, we provide Schur() also for traditional matrices: a.m <- function(M) unname(as(M, "matrix")) a <- a.m(A) Sch.a <- Schur(a) stopifnot(identical(Sch.a, list(Q = a.m(Sch.A @ Q), T = a.m(Sch.A @ T), EValues = Sch.A@EValues)), all.equal(a, with(Sch.a, Q %*% T %*% t(Q))) )