Title: | Complex Matrix Algebra with 'Eigen' |
---|---|
Description: | Matrix algebra using the 'Eigen' C++ library: determinant, rank, inverse, pseudo-inverse, kernel and image, QR decomposition, Cholesky decomposition, Schur decomposition, Hessenberg decomposition, linear least-squares problems. Also provides matrix functions such as exponential, logarithm, power, sine and cosine. Complex matrices are supported. |
Authors: | Stéphane Laurent |
Maintainer: | Stéphane Laurent <[email protected]> |
License: | GPL-3 |
Version: | 1.3.0 |
Built: | 2024-10-25 04:46:20 UTC |
Source: | https://github.com/stla/eigenr |
Absolute value of the determinant of a real matrix.
Eigen_absdet(M)
Eigen_absdet(M)
M |
a real square matrix |
The absolute value of the determinant of M
.
'Eigen_absdet(M)' is not faster than 'abs(Eigen_det(M))'.
set.seed(666L) M <- matrix(rpois(25L, 1), 5L, 5L) Eigen_absdet(M)
set.seed(666L) M <- matrix(rpois(25L, 1), 5L, 5L) Eigen_absdet(M)
Cholesky decomposition of a symmetric or Hermitian matrix.
Eigen_chol(M)
Eigen_chol(M)
M |
a square symmetric/Hermitian positive-definite matrix or
|
Symmetry is not checked; only the lower triangular part of
M
is used.
The upper triangular factor of the Cholesky decomposition of
M
.
M <- rbind(c(5,1), c(1,3)) U <- Eigen_chol(M) t(U) %*% U # this is `M` # a Hermitian example: A <- rbind(c(1,1i), c(1i,2)) ( M <- A %*% t(Conj(A)) ) try(chol(M)) # fails U <- Eigen_chol(M) t(Conj(U)) %*% U # this is `M` # a sparse example M <- asSparseMatrix(diag(1:5)) Eigen_chol(M)
M <- rbind(c(5,1), c(1,3)) U <- Eigen_chol(M) t(U) %*% U # this is `M` # a Hermitian example: A <- rbind(c(1,1i), c(1i,2)) ( M <- A %*% t(Conj(A)) ) try(chol(M)) # fails U <- Eigen_chol(M) t(Conj(U)) %*% U # this is `M` # a sparse example M <- asSparseMatrix(diag(1:5)) Eigen_chol(M)
Complex Schur decomposition of a square matrix.
Eigen_complexSchur(M)
Eigen_complexSchur(M)
M |
real or complex square matrix |
See Eigen::ComplexSchur.
A list with the T
and U
matrices.
library(EigenR) M <- cbind(c(3, 2i, 1+3i), c(1, 1i, 1), c(5, 0, -2i)) schur <- Eigen_complexSchur(M) T <- schur$T U <- schur$U M - U %*% T %*% t(Conj(U))
library(EigenR) M <- cbind(c(3, 2i, 1+3i), c(1, 1i, 1), c(5, 0, -2i)) schur <- Eigen_complexSchur(M) T <- schur$T U <- schur$U M - U %*% T %*% t(Conj(U))
Matrix cosine of a real or complex square matrix.
Eigen_cos(M)
Eigen_cos(M)
M |
a square matrix, real or complex |
The matrix cosine of M
.
library(EigenR) M <- toeplitz(c(1,2,3)) cosM <- Eigen_cos(M) sinM <- Eigen_sin(M) cosM %*% cosM + sinM %*% sinM # identity matrix
library(EigenR) M <- toeplitz(c(1,2,3)) cosM <- Eigen_cos(M) sinM <- Eigen_sin(M) cosM %*% cosM + sinM %*% sinM # identity matrix
Matrix hyperbolic cosine of a real or complex square matrix.
Eigen_cosh(M)
Eigen_cosh(M)
M |
a square matrix, real or complex |
The matrix hyperbolic cosine of M
.
library(EigenR) M <- toeplitz(c(1,2,3)) Eigen_cosh(M) (Eigen_exp(M) + Eigen_exp(-M)) / 2 # identical
library(EigenR) M <- toeplitz(c(1,2,3)) Eigen_cosh(M) (Eigen_exp(M) + Eigen_exp(-M)) / 2 # identical
Determinant of a real or complex matrix.
Eigen_det(M)
Eigen_det(M)
M |
a square matrix or |
The determinant of M
.
set.seed(666) M <- matrix(rpois(25, 1), 5L, 5L) Eigen_det(M) # determinants of complex matrices are supported: Eigen_det(M + 1i * M) # as well as determinants of sparse matrices: Eigen_det(asSparseMatrix(M)) Eigen_det(asSparseMatrix(M + 1i * M))
set.seed(666) M <- matrix(rpois(25, 1), 5L, 5L) Eigen_det(M) # determinants of complex matrices are supported: Eigen_det(M + 1i * M) # as well as determinants of sparse matrices: Eigen_det(asSparseMatrix(M)) Eigen_det(asSparseMatrix(M + 1i * M))
Exponential of a real or complex square matrix.
Eigen_exp(M)
Eigen_exp(M)
M |
a square matrix, real or complex |
The exponential of M
.
Hessenberg decomposition of a square matrix.
Eigen_Hessenberg(M)
Eigen_Hessenberg(M)
M |
real or complex square matrix |
See Eigen::HessenbergDecomposition.
A list with the H
and Q
matrices.
library(EigenR) M <- cbind(c(3, 2i, 1+3i), c(1, 1i, 1), c(5, 0, -2i)) Eigen_Hessenberg(M)
library(EigenR) M <- cbind(c(3, 2i, 1+3i), c(1, 1i, 1), c(5, 0, -2i)) Eigen_Hessenberg(M)
Inverse of a real or complex matrix.
Eigen_inverse(M)
Eigen_inverse(M)
M |
an invertible square matrix, real or complex |
The inverse matrix of M
.
Checks whether a matrix represents an injective linear map (i.e. has trivial kernel).
Eigen_isInjective(M)
Eigen_isInjective(M)
M |
a matrix, real or complex |
A Boolean value indicating whether M
represents an injective
linear map.
set.seed(666L) M <- matrix(rpois(35L, 1), 5L, 7L) Eigen_isInjective(M)
set.seed(666L) M <- matrix(rpois(35L, 1), 5L, 7L) Eigen_isInjective(M)
Checks whether a matrix is invertible.
Eigen_isInvertible(M)
Eigen_isInvertible(M)
M |
a matrix, real or complex |
A Boolean value indicating whether M
is invertible.
set.seed(666L) M <- matrix(rpois(25L, 1), 5L, 5L) Eigen_isInvertible(M)
set.seed(666L) M <- matrix(rpois(25L, 1), 5L, 5L) Eigen_isInvertible(M)
Checks whether a matrix represents a surjective linear map.
Eigen_isSurjective(M)
Eigen_isSurjective(M)
M |
a matrix, real or complex |
A Boolean value indicating whether M
represents a surjective
linear map.
set.seed(666L) M <- matrix(rpois(35L, 1), 7L, 5L) Eigen_isSurjective(M)
set.seed(666L) M <- matrix(rpois(35L, 1), 7L, 5L) Eigen_isSurjective(M)
Kernel (null-space) of a real or complex matrix.
Eigen_kernel(M, method = "COD")
Eigen_kernel(M, method = "COD")
M |
a matrix, real or complex |
method |
one of |
A basis of the kernel of M
. With method = "COD"
, the
basis is orthonormal, while it is not with method = "LU"
.
set.seed(666) M <- matrix(rgamma(30L, 12, 1), 10L, 3L) M <- cbind(M, M[,1]+M[,2], M[,2]+2*M[,3]) # basis of the kernel of `M`: Eigen_kernel(M, method = "LU") # orthonormal basis of the kernel of `M`: Eigen_kernel(M, method = "COD")
set.seed(666) M <- matrix(rgamma(30L, 12, 1), 10L, 3L) M <- cbind(M, M[,1]+M[,2], M[,2]+2*M[,3]) # basis of the kernel of `M`: Eigen_kernel(M, method = "LU") # orthonormal basis of the kernel of `M`: Eigen_kernel(M, method = "COD")
Dimension of the kernel of a matrix.
Eigen_kernelDimension(M)
Eigen_kernelDimension(M)
M |
a matrix, real or complex |
An integer, the dimension of the kernel of M
.
Eigen_isInjective
, Eigen_kernel
.
set.seed(666L) M <- matrix(rpois(35L, 1), 5L, 7L) Eigen_kernelDimension(M)
set.seed(666L) M <- matrix(rpois(35L, 1), 5L, 7L) Eigen_kernelDimension(M)
Logarithm of a real or complex square matrix, when possible.
Eigen_log(M)
Eigen_log(M)
M |
a square matrix, real or complex |
The logarithm of a matrix does not always exist. See matrix logarithm.
The logarithm of M
.
Logarithm of the absolute value of the determinant of a real matrix.
Eigen_logabsdet(M)
Eigen_logabsdet(M)
M |
a real square matrix |
The logarithm of the absolute value of the determinant of M
.
'Eigen_logabsdet(M)' is not faster than 'log(abs(Eigen_det(M)))'.
set.seed(666L) M <- matrix(rpois(25L, 1), 5L, 5L) Eigen_logabsdet(M)
set.seed(666L) M <- matrix(rpois(25L, 1), 5L, 5L) Eigen_logabsdet(M)
Solves a linear least-squares problem.
Eigen_lsSolve(A, b, method = "cod")
Eigen_lsSolve(A, b, method = "cod")
A |
a |
b |
a vector of length |
method |
the method used to solve the problem, either |
The solution X
of the least-squares problem AX ~= b
(similar to lm.fit(A, b)$coefficients
). This is a matrix if
b
is a matrix, or a vector if b
is a vector.
set.seed(129) n <- 7; p <- 2 A <- matrix(rnorm(n * p), n, p) b <- rnorm(n) lsfit <- Eigen_lsSolve(A, b) b - A %*% lsfit # residuals
set.seed(129) n <- 7; p <- 2 A <- matrix(rnorm(n * p), n, p) b <- rnorm(n) lsfit <- Eigen_lsSolve(A, b) b - A %*% lsfit # residuals
Pseudo-inverse of a real or complex matrix (Moore-Penrose generalized inverse).
Eigen_pinverse(M)
Eigen_pinverse(M)
M |
a matrix, real or complex, not necessarily square |
The pseudo-inverse matrix of M
.
library(EigenR) M <- rbind( toeplitz(c(3, 2, 1)), toeplitz(c(4, 5, 6)) ) Mplus <- Eigen_pinverse(M) all.equal(M, M %*% Mplus %*% M) all.equal(Mplus, Mplus %*% M %*% Mplus) #' a complex matrix A <- M + 1i * M[, c(3L, 2L, 1L)] Aplus <- Eigen_pinverse(A) AAplus <- A %*% Aplus all.equal(AAplus, t(Conj(AAplus))) #' `A %*% Aplus` is Hermitian AplusA <- Aplus %*% A all.equal(AplusA, t(Conj(AplusA))) #' `Aplus %*% A` is Hermitian
library(EigenR) M <- rbind( toeplitz(c(3, 2, 1)), toeplitz(c(4, 5, 6)) ) Mplus <- Eigen_pinverse(M) all.equal(M, M %*% Mplus %*% M) all.equal(Mplus, Mplus %*% M %*% Mplus) #' a complex matrix A <- M + 1i * M[, c(3L, 2L, 1L)] Aplus <- Eigen_pinverse(A) AAplus <- A %*% Aplus all.equal(AAplus, t(Conj(AAplus))) #' `A %*% Aplus` is Hermitian AplusA <- Aplus %*% A all.equal(AplusA, t(Conj(AplusA))) #' `Aplus %*% A` is Hermitian
Matricial power of a real or complex square matrix, when possible.
Eigen_pow(M, p)
Eigen_pow(M, p)
M |
a square matrix, real or complex |
p |
a number, real or complex, the power exponent |
The power is defined with the help of the exponential and the logarithm. See matrix power.
The matrix M
raised at the power p
.
QR decomposition of a real or complex matrix.
Eigen_QR(M)
Eigen_QR(M)
M |
a matrix, real or complex |
A list with the Q
matrix and the R
matrix.
M <- cbind(c(1,2,3), c(4,5,6)) x <- Eigen_QR(M) x$Q %*% x$R
M <- cbind(c(1,2,3), c(4,5,6)) x <- Eigen_QR(M) x$Q %*% x$R
Range (column-space, image, span) of a real or complex matrix.
Eigen_range(M, method = "QR")
Eigen_range(M, method = "QR")
M |
a matrix, real or complex |
method |
one of |
A basis of the range of M
. With method = "LU"
, the
basis is not orthonormal, while it is with method = "QR"
and
method = "COD"
.
Rank of a real or complex matrix.
Eigen_rank(M)
Eigen_rank(M)
M |
a matrix, real or complex |
The rank of M
.
Real Schur decomposition of a square matrix.
Eigen_realSchur(M)
Eigen_realSchur(M)
M |
real square matrix |
See Eigen::RealSchur.
A list with the T
and U
matrices.
library(EigenR) M <- cbind(c(3, 2, 3), c(1, 1, 1), c(5, 0, -2)) schur <- Eigen_realSchur(M) T <- schur$T U <- schur$U M - U %*% T %*% t(U)
library(EigenR) M <- cbind(c(3, 2, 3), c(1, 1, 1), c(5, 0, -2)) schur <- Eigen_realSchur(M) T <- schur$T U <- schur$U M - U %*% T %*% t(U)
Matrix sine of a real or complex square matrix.
Eigen_sin(M)
Eigen_sin(M)
M |
a square matrix, real or complex |
The matrix sine of M
.
Matrix hyperbolic sine of a real or complex square matrix.
Eigen_sinh(M)
Eigen_sinh(M)
M |
a square matrix, real or complex |
The matrix hyperbolic sine of M
.
library(EigenR) M <- toeplitz(c(1,2,3)) Eigen_sinh(M) (Eigen_exp(M) - Eigen_exp(-M)) / 2 # identical
library(EigenR) M <- toeplitz(c(1,2,3)) Eigen_sinh(M) (Eigen_exp(M) - Eigen_exp(-M)) / 2 # identical
Square root of a real or complex square matrix, when possible.
Eigen_sqrt(M)
Eigen_sqrt(M)
M |
a square matrix, real or complex |
See matrix square root.
A square root of M
.
# Rotation matrix over 60 degrees: M <- cbind(c(cos(pi/3), sin(pi/3)), c(-sin(pi/3), cos(pi/3))) # Its square root, the rotation matrix over 30 degrees: Eigen_sqrt(M)
# Rotation matrix over 60 degrees: M <- cbind(c(cos(pi/3), sin(pi/3)), c(-sin(pi/3), cos(pi/3))) # Its square root, the rotation matrix over 30 degrees: Eigen_sqrt(M)
Cholesky-'UtDU' decomposition of a symmetric or Hermitian matrix.
Eigen_UtDU(M)
Eigen_UtDU(M)
M |
a square symmetric/Hermitian positive or negative semidefinite matrix, real/complex |
Symmetry is not checked; only the lower triangular part of
M
is used.
The Cholesky-'UtDU' decomposition of M
in a list
(see example).
x <- matrix(c(1:5, (1:5)^2), 5, 2) x <- cbind(x, x[, 1] + 3*x[, 2]) M <- crossprod(x) UtDU <- Eigen_UtDU(M) U <- UtDU$U D <- UtDU$D perm <- UtDU$perm UP <- U[, perm] t(UP) %*% diag(D) %*% UP # this is `M`
x <- matrix(c(1:5, (1:5)^2), 5, 2) x <- cbind(x, x[, 1] + 3*x[, 2]) M <- crossprod(x) UtDU <- Eigen_UtDU(M) U <- UtDU$U D <- UtDU$D perm <- UtDU$perm UP <- U[, perm] t(UP) %*% diag(D) %*% UP # this is `M`
Constructs a sparse matrix, real or complex.
SparseMatrix(i, j, Mij, nrows, ncols) ## S3 method for class 'SparseMatrix' print(x, ...) asSparseMatrix(M)
SparseMatrix(i, j, Mij, nrows, ncols) ## S3 method for class 'SparseMatrix' print(x, ...) asSparseMatrix(M)
i , j
|
indices of the non-zero coefficients |
Mij |
values of the non-zero coefficients; must be a vector of the same
length as |
nrows , ncols
|
dimensions of the matrix |
x |
a |
... |
ignored |
M |
a matrix, real or complex |
A list with the class SparseMatrix
.
set.seed(666) ( M <- matrix(rpois(50L, 1), 10L, 5L) ) asSparseMatrix(M)
set.seed(666) ( M <- matrix(rpois(50L, 1), 10L, 5L) ) asSparseMatrix(M)