Title: | Tools for Tensor Analysis and Decomposition |
---|---|
Description: | A set of tools for creation, manipulation, and modeling of tensors with arbitrary number of modes. A tensor in the context of data analysis is a multidimensional array. rTensor does this by providing a S4 class 'Tensor' that wraps around the base 'array' class. rTensor provides common tensor operations as methods, including matrix unfolding, summing/averaging across modes, calculating the Frobenius norm, and taking the inner product between two tensors. Familiar array operations are overloaded, such as index subsetting via '[' and element-wise operations. rTensor also implements various tensor decomposition, including CP, GLRAM, MPCA, PVD, and Tucker. For tensors with 3 modes, rTensor also implements transpose, t-product, and t-SVD, as defined in Kilmer et al. (2013). Some auxiliary functions include the Khatri-Rao product, Kronecker product, and the Hadamard product for a list of matrices. |
Authors: | James Li and Jacob Bien and Martin Wells |
Maintainer: | Koki Tsuyuzaki <[email protected]> |
License: | MIT + file LICENSE |
Version: | 1.4.8 |
Built: | 2024-11-01 11:15:24 UTC |
Source: | https://github.com/rikenbit/rtensor |
This package is centered around the Tensor-class
, which defines a S4 class for tensors of arbitrary number of modes. A vignette and/or a possible paper will be included in a future release of this package.
This page will summarize the full functionality of this package. Note that since all the methods associated with S4 class Tensor-class
are documented there, we will not duplicate it here.
The remaining functions can be split into two groups: the first is a set of tensor decompositions, and the second is a set of helper functions that are useful in tensor manipulation.
rTensor implements the following tensor decompositions:
cp
Canonical Polyadic (CP) decomposition
tucker
General Tucker decomposition
mpca
Multilinear Principal Component Analysis; note that for 3-Tensors this is also known as Generalized Low Rank Approximation of Matrices(GLRAM)
hosvd
(Truncated-)Higher-order singular value decomposition
t_svd
Tensor singular value decomposition; 3-Tensors only; also note that there is an asociated reconstruction function t_svd_reconstruct
pvd
Population value decomposition of images; 3-Tensors only
rTensor also provides a set functions for tensors multiplication:
ttm
Tensor times matrix, aka m-mode product
ttl
Tensor times list (of matrices)
t_mult
Tensor product based on block circulant unfolding; only implemented for a pair of 3-Tensors
...as well as for matrices:
hadamard_list
Computes the Hadamard (element-wise) product of a list of matrices
kronecker_list
Computes the Kronecker product of a list of matrices
khatri_rao
Computes the Khatri-Rao product of two matrices
khatri_rao_list
Computes the Khatri-Rao product of a list of matrices
fold
General folding of a matrix into a tensor
k_fold
Inverse operation for k_unfold
unmatvec
Inverse operation for matvec
For more information on any of the functions, please consult the individual man pages.
James Li [email protected], Jacob Bien, and Martin T. Wells
Extends '[' and '[<-' from the base array class for the Tensor class. Works exactly as it would for the base 'array' class.
## S4 method for signature 'Tensor' x[i, j, ..., drop = TRUE] ## S4 replacement method for signature 'Tensor' x[i, j, ...] <- value
## S4 method for signature 'Tensor' x[i, j, ..., drop = TRUE] ## S4 replacement method for signature 'Tensor' x[i, j, ...] <- value
x |
Tensor to be subset |
i , j , ...
|
indices that specify the extents of the sub-tensor |
drop |
whether or not to reduce the number of modes to exclude those that have '1' as the mode |
value |
either vector, matrix, or array that will replace the subtensor |
x[i,j,...,drop=TRUE]
an object of class Tensor
tnsr <- rand_tensor() tnsr[1,2,3] tnsr[3,1,] tnsr[,,5] tnsr[,,5,drop=FALSE] tnsr[1,2,3] <- 3; tnsr[1,2,3] tnsr[3,1,] <- rep(0,5); tnsr[3,1,] tnsr[,2,] <- matrix(0,nrow=3,ncol=5); tnsr[,2,]
tnsr <- rand_tensor() tnsr[1,2,3] tnsr[3,1,] tnsr[,,5] tnsr[,,5,drop=FALSE] tnsr[1,2,3] <- 3; tnsr[1,2,3] tnsr[3,1,] <- rep(0,5); tnsr[3,1,] tnsr[,2,] <- matrix(0,nrow=3,ncol=5); tnsr[,2,]
Create a Tensor-class
object from an array
, matrix
, or vector
.
as.tensor(x, drop = FALSE)
as.tensor(x, drop = FALSE)
x |
an instance of |
drop |
whether or not modes of 1 should be dropped |
a Tensor-class
object
#From vector vec <- runif(3); vecT <- as.tensor(vec); vecT #From matrix mat <- matrix(runif(2*3),nrow=2,ncol=3) matT <- as.tensor(mat); matT #From array indices <- c(2,3,4) arr <- array(runif(prod(indices)), dim = indices) arrT <- as.tensor(arr); arrT
#From vector vec <- runif(3); vecT <- as.tensor(vec); vecT #From matrix mat <- matrix(runif(2*3),nrow=2,ncol=3) matT <- as.tensor(mat); matT #From array indices <- c(2,3,4) arr <- array(runif(prod(indices)), dim = indices) arrT <- as.tensor(arr); arrT
Canonical Polyadic (CP) decomposition of a tensor, aka CANDECOMP/PARAFRAC. Approximate a K-Tensor using a sum of num_components
rank-1 K-Tensors. A rank-1 K-Tensor can be written as an outer product of K vectors. There are a total of num_compoents *tnsr@num_modes
vectors in the output, stored in tnsr@num_modes
matrices, each with num_components
columns. This is an iterative algorithm, with two possible stopping conditions: either relative error in Frobenius norm has gotten below tol
, or the max_iter
number of iterations has been reached. For more details on CP decomposition, consult Kolda and Bader (2009).
cp(tnsr, num_components = NULL, max_iter = 25, tol = 1e-05)
cp(tnsr, num_components = NULL, max_iter = 25, tol = 1e-05)
tnsr |
Tensor with K modes |
num_components |
the number of rank-1 K-Tensors to use in approximation |
max_iter |
maximum number of iterations if error stays above |
tol |
relative Frobenius norm error tolerance |
Uses the Alternating Least Squares (ALS) estimation procedure. A progress bar is included to help monitor operations on large tensors.
a list containing the following
lambdas
a vector of normalizing constants, one for each component
U
a list of matrices - one for each mode - each matrix with num_components
columns
conv
whether or not resid
< tol
by the last iteration
norm_percent
the percent of Frobenius norm explained by the approximation
est
estimate of tnsr
after compression
fnorm_resid
the Frobenius norm of the error fnorm(est-tnsr)
all_resids
vector containing the Frobenius norm of error for all the iterations
T. Kolda, B. Bader, "Tensor decomposition and applications". SIAM Applied Mathematics and Applications 2009.
### How to retrieve faces_tnsr from figshare # faces_tnsr <- load_orl() # subject <- faces_tnsr[,,14,] dummy_faces_tnsr <- rand_tensor(c(92,112,40,10)) subject <- dummy_faces_tnsr[,,14,] cpD <- cp(subject, num_components=3) cpD$conv cpD$norm_percent plot(cpD$all_resids)
### How to retrieve faces_tnsr from figshare # faces_tnsr <- load_orl() # subject <- faces_tnsr[,,14,] dummy_faces_tnsr <- rand_tensor(c(92,112,40,10)) subject <- dummy_faces_tnsr[,,14,] cpD <- cp(subject, num_components=3) cpD$conv cpD$norm_percent plot(cpD$all_resids)
DEPRECATED. Please see unmatvec
cs_fold(mat, m = NULL, modes = NULL)
cs_fold(mat, m = NULL, modes = NULL)
mat |
matrix to be folded |
m |
the mode corresponding to cs_unfold |
modes |
the original modes of the tensor |
DEPRECATED. Please see matvec-methods
and unfold-methods
.
cs_unfold(tnsr, m) ## S4 method for signature 'Tensor' cs_unfold(tnsr, m = NULL)
cs_unfold(tnsr, m) ## S4 method for signature 'Tensor' cs_unfold(tnsr, m = NULL)
tnsr |
Tensor instance |
m |
mode to be unfolded on |
cs_unfold(tnsr,m=NULL)
Return the vector of modes from a tensor
## S4 method for signature 'Tensor' dim(x)
## S4 method for signature 'Tensor' dim(x)
x |
the Tensor instance |
dim(x)
an integer vector of the modes associated with x
tnsr <- rand_tensor() dim(tnsr)
tnsr <- rand_tensor() dim(tnsr)
Returns the Frobenius norm of the Tensor instance.
fnorm(tnsr) ## S4 method for signature 'Tensor' fnorm(tnsr)
fnorm(tnsr) ## S4 method for signature 'Tensor' fnorm(tnsr)
tnsr |
the Tensor instance |
fnorm(tnsr)
numeric Frobenius norm of x
tnsr <- rand_tensor() fnorm(tnsr)
tnsr <- rand_tensor() fnorm(tnsr)
General folding of a matrix into a Tensor. This is designed to be the inverse function to unfold-methods
, with the same ordering of the indices. This amounts to following: if we were to unfold a Tensor using a set of row_idx
and col_idx
, then we can fold the resulting matrix back into the original Tensor using the same row_idx
and col_idx
.
fold(mat, row_idx = NULL, col_idx = NULL, modes = NULL)
fold(mat, row_idx = NULL, col_idx = NULL, modes = NULL)
mat |
matrix to be folded into a Tensor |
row_idx |
the indices of the modes that are mapped onto the row space |
col_idx |
the indices of the modes that are mapped onto the column space |
modes |
the modes of the output Tensor |
This function uses aperm
as the primary workhorse.
Tensor object with modes given by modes
T. Kolda, B. Bader, "Tensor decomposition and applications". SIAM Applied Mathematics and Applications 2009.
unfold-methods
, k_fold
, unmatvec
tnsr <- new("Tensor",3L,c(3L,4L,5L),data=runif(60)) matT3<-unfold(tnsr,row_idx=2,col_idx=c(3,1)) identical(fold(matT3,row_idx=2,col_idx=c(3,1),modes=c(3,4,5)),tnsr)
tnsr <- new("Tensor",3L,c(3L,4L,5L),data=runif(60)) matT3<-unfold(tnsr,row_idx=2,col_idx=c(3,1)) identical(fold(matT3,row_idx=2,col_idx=c(3,1),modes=c(3,4,5)),tnsr)
Returns the hadamard (element-wise) product from a list of matrices or vectors. Commonly used for n-mode products and various Tensor decompositions.
hadamard_list(L)
hadamard_list(L)
L |
list of matrices or vectors |
matrix that is the hadamard product
The modes/dimensions of each element in the list must match.
kronecker_list
, khatri_rao_list
lizt <- list('mat1' = matrix(runif(40),ncol=4), 'mat2' = matrix(runif(40),ncol=4), 'mat3' = matrix(runif(40),ncol=4)) dim(hadamard_list(lizt))
lizt <- list('mat1' = matrix(runif(40),ncol=4), 'mat2' = matrix(runif(40),ncol=4), 'mat3' = matrix(runif(40),ncol=4)) dim(hadamard_list(lizt))
Extend head for Tensor
## S4 method for signature 'Tensor' head(x, ...)
## S4 method for signature 'Tensor' head(x, ...)
x |
the Tensor instance |
... |
additional parameters to be passed into head() |
head(x,...)
tnsr <- rand_tensor() head(tnsr)
tnsr <- rand_tensor() head(tnsr)
Higher-order SVD of a K-Tensor. Write the K-Tensor as a (m-mode) product of a core Tensor (possibly smaller modes) and K orthogonal factor matrices. Truncations can be specified via ranks
(making them smaller than the original modes of the K-Tensor will result in a truncation). For the mathematical details on HOSVD, consult Lathauwer et. al. (2000).
hosvd(tnsr, ranks = NULL)
hosvd(tnsr, ranks = NULL)
tnsr |
Tensor with K modes |
ranks |
a vector of desired modes in the output core tensor, default is |
A progress bar is included to help monitor operations on large tensors.
a list containing the following:
Z
core tensor with modes speficied by ranks
U
a list of orthogonal matrices, one for each mode
est
estimate of tnsr
after compression
fnorm_resid
the Frobenius norm of the error fnorm(est-tnsr)
- if there was no truncation, then this is on the order of mach_eps * fnorm.
The length of ranks
must match tnsr@num_modes
.
L. Lathauwer, B.Moor, J. Vanderwalle "A multilinear singular value decomposition". Journal of Matrix Analysis and Applications 2000.
tnsr <- rand_tensor(c(6,7,8)) hosvdD <- hosvd(tnsr) plot(hosvdD$fnorm_resid) hosvdD2 <- hosvd(tnsr,ranks=c(3,3,4)) plot(hosvdD2$fnorm_resid)
tnsr <- rand_tensor(c(6,7,8)) hosvdD <- hosvd(tnsr) plot(hosvdD$fnorm_resid) hosvdD2 <- hosvd(tnsr,ranks=c(3,3,4)) plot(hosvdD2$fnorm_resid)
Not designed to be called by the user. Use as.tensor
instead.
## S4 method for signature 'Tensor' initialize(.Object, num_modes = NULL, modes = NULL, data = NULL)
## S4 method for signature 'Tensor' initialize(.Object, num_modes = NULL, modes = NULL, data = NULL)
.Object |
the tensor object |
num_modes |
number of modes of the tensor |
modes |
modes of the tensor |
data |
can be vector, matrix, or array |
as.tensor
Returns the inner product between two Tensors
innerProd(tnsr1, tnsr2) ## S4 method for signature 'Tensor,Tensor' innerProd(tnsr1, tnsr2)
innerProd(tnsr1, tnsr2) ## S4 method for signature 'Tensor,Tensor' innerProd(tnsr1, tnsr2)
tnsr1 |
first Tensor instance |
tnsr2 |
second Tensor instance |
innerProd(tnsr1,tnsr2)
inner product between x1
and x2
tnsr1 <- rand_tensor() tnsr2 <- rand_tensor() innerProd(tnsr1,tnsr2)
tnsr1 <- rand_tensor() tnsr2 <- rand_tensor() innerProd(tnsr1,tnsr2)
k-mode folding of a matrix into a Tensor. This is the inverse funtion to k_unfold
in the m mode. In particular, k_fold(k_unfold(tnsr, m),m,getModes(tnsr))
will result in the original Tensor.
k_fold(mat, m = NULL, modes = NULL)
k_fold(mat, m = NULL, modes = NULL)
mat |
matrix to be folded into a Tensor |
m |
the index of the mode that is mapped onto the row indices |
modes |
the modes of the output Tensor |
This is a wrapper function to fold
.
Tensor object with modes given by modes
T. Kolda, B. Bader, "Tensor decomposition and applications". SIAM Applied Mathematics and Applications 2009.
k_unfold-methods
, fold
, unmatvec
tnsr <- new("Tensor",3L,c(3L,4L,5L),data=runif(60)) matT2<-k_unfold(tnsr,m=2) identical(k_fold(matT2,m=2,modes=c(3,4,5)),tnsr)
tnsr <- new("Tensor",3L,c(3L,4L,5L),data=runif(60)) matT2<-k_unfold(tnsr,m=2) identical(k_fold(matT2,m=2,modes=c(3,4,5)),tnsr)
Unfolding of a tensor by mapping the kth mode (specified through parameter m
), and all other modes onto the column space. This the most common type of unfolding operation for Tucker decompositions and its variants. Also known as k-mode matricization.
k_unfold(tnsr, m) ## S4 method for signature 'Tensor' k_unfold(tnsr, m = NULL)
k_unfold(tnsr, m) ## S4 method for signature 'Tensor' k_unfold(tnsr, m = NULL)
tnsr |
the Tensor instance |
m |
the index of the mode to unfold on |
k_unfold(tnsr,m=NULL)
matrix with x@modes[m]
rows and prod(x@modes[-m])
columns
T. Kolda and B. Bader, "Tensor decomposition and applications". SIAM Applied Mathematics and Applications 2009.
matvec-methods
and unfold-methods
tnsr <- rand_tensor() matT2<-rs_unfold(tnsr,m=2)
tnsr <- rand_tensor() matT2<-rs_unfold(tnsr,m=2)
Returns the Khatri-Rao (column-wise Kronecker) product of two matrices. If the inputs are vectors then this is the same as the Kronecker product.
khatri_rao(x, y)
khatri_rao(x, y)
x |
first matrix |
y |
second matrix |
matrix that is the Khatri-Rao product
The number of columns must match in the two inputs.
dim(khatri_rao(matrix(runif(12),ncol=4),matrix(runif(12),ncol=4)))
dim(khatri_rao(matrix(runif(12),ncol=4),matrix(runif(12),ncol=4)))
Returns the Khatri-Rao product from a list of matrices or vectors. Commonly used for n-mode products and various Tensor decompositions.
khatri_rao_list(L, reverse = FALSE)
khatri_rao_list(L, reverse = FALSE)
L |
list of matrices or vectors |
reverse |
whether or not to reverse the order |
matrix that is the Khatri-Rao product
The number of columns must match in every element of the input list.
smalllizt <- list('mat1' = matrix(runif(12),ncol=4), 'mat2' = matrix(runif(12),ncol=4), 'mat3' = matrix(runif(12),ncol=4)) dim(khatri_rao_list(smalllizt))
smalllizt <- list('mat1' = matrix(runif(12),ncol=4), 'mat2' = matrix(runif(12),ncol=4), 'mat3' = matrix(runif(12),ncol=4)) dim(khatri_rao_list(smalllizt))
Returns the Kronecker product from a list of matrices or vectors. Commonly used for n-mode products and various Tensor decompositions.
kronecker_list(L)
kronecker_list(L)
L |
list of matrices or vectors |
matrix that is the Kronecker product
hadamard_list
, khatri_rao_list
, kronecker
smalllizt <- list('mat1' = matrix(runif(12),ncol=4), 'mat2' = matrix(runif(12),ncol=4), 'mat3' = matrix(runif(12),ncol=4)) dim(kronecker_list(smalllizt))
smalllizt <- list('mat1' = matrix(runif(12),ncol=4), 'mat2' = matrix(runif(12),ncol=4), 'mat3' = matrix(runif(12),ncol=4)) dim(kronecker_list(smalllizt))
A dataset containing pictures of 40 individuals under 10 different lightings. Each image has 92 x 112 pixels. Structured as a 4-tensor with modes 92 x 112 x 40 x 10. The data is now stored in figshare https://ndownloader.figshare.com/files/28005669
load_orl()
load_orl()
A Tensor object with modes 92 x 112 x 40 x 10. The first two modes correspond to the image pixels, the third mode corresponds to the individual, and the last mode correpsonds to the lighting.
https://www.kaggle.com/kasikrit/att-database-of-faces
AT&T Laboratories Cambridge. https://www.kaggle.com/kasikrit/att-database-of-faces
F. Samaria, A. Harter, "Parameterisation of a Stochastic Model for Human Face Identification". IEEE Workshop on Applications of Computer Vision 1994.
For 3-tensors only. Stacks the slices along the third mode. This is the prevalent unfolding for T-SVD and T-MULT based on block circulant matrices.
matvec(tnsr) ## S4 method for signature 'Tensor' matvec(tnsr)
matvec(tnsr) ## S4 method for signature 'Tensor' matvec(tnsr)
tnsr |
the Tensor instance |
matvec(tnsr)
matrix with prod(x@modes[-m])
rows and x@modes[m]
columns
M. Kilmer, K. Braman, N. Hao, and R. Hoover, "Third-order tensors as operators on matrices: a theoretical and computational framework with applications in imaging". SIAM Journal on Matrix Analysis and Applications 2013.
k_unfold-methods
and unfold-methods
tnsr <- rand_tensor(c(2,3,4)) matT1<- matvec(tnsr)
tnsr <- rand_tensor(c(2,3,4)) matT1<- matvec(tnsr)
Given a mode for a K-tensor, this returns the K-1 tensor resulting from taking the mean across that particular mode.
modeMean(tnsr, m, drop) ## S4 method for signature 'Tensor' modeMean(tnsr, m = NULL, drop = FALSE)
modeMean(tnsr, m, drop) ## S4 method for signature 'Tensor' modeMean(tnsr, m = NULL, drop = FALSE)
tnsr |
the Tensor instance |
m |
the index of the mode to average across |
drop |
whether or not mode m should be dropped |
modeMean(tnsr,m=NULL,drop=FALSE)
K-1 or K Tensor, where K = x@num_modes
tnsr <- rand_tensor() modeMean(tnsr,1,drop=TRUE)
tnsr <- rand_tensor() modeMean(tnsr,1,drop=TRUE)
Given a mode for a K-tensor, this returns the K-1 tensor resulting from summing across that particular mode.
modeSum(tnsr, m, drop) ## S4 method for signature 'Tensor' modeSum(tnsr, m = NULL, drop = FALSE)
modeSum(tnsr, m, drop) ## S4 method for signature 'Tensor' modeSum(tnsr, m = NULL, drop = FALSE)
tnsr |
the Tensor instance |
m |
the index of the mode to sum across |
drop |
whether or not mode m should be dropped |
modeSum(tnsr,m=NULL,drop=FALSE)
K-1 or K tensor, where K = x@num_modes
tnsr <- rand_tensor() modeSum(tnsr,3,drop=TRUE)
tnsr <- rand_tensor() modeSum(tnsr,3,drop=TRUE)
This is basically the Tucker decomposition of a K-Tensor, tucker
, with one of the modes uncompressed. If K = 3, then this is also known as the Generalized Low Rank Approximation of Matrices (GLRAM). This implementation assumes that the last mode is the measurement mode and hence uncompressed. This is an iterative algorithm, with two possible stopping conditions: either relative error in Frobenius norm has gotten below tol
, or the max_iter
number of iterations has been reached. For more details on the MPCA of tensors, consult Lu et al. (2008).
mpca(tnsr, ranks = NULL, max_iter = 25, tol = 1e-05)
mpca(tnsr, ranks = NULL, max_iter = 25, tol = 1e-05)
tnsr |
Tensor with K modes |
ranks |
a vector of the compressed modes of the output core Tensor, this has length K-1 |
max_iter |
maximum number of iterations if error stays above |
tol |
relative Frobenius norm error tolerance |
Uses the Alternating Least Squares (ALS) estimation procedure. A progress bar is included to help monitor operations on large tensors.
a list containing the following:
Z_ext
the extended core tensor, with the first K-1 modes given by ranks
U
a list of K-1 orthgonal factor matrices - one for each compressed mode, with the number of columns of the matrices given by ranks
conv
whether or not resid
< tol
by the last iteration
est
estimate of tnsr
after compression
norm_percent
the percent of Frobenius norm explained by the approximation
fnorm_resid
the Frobenius norm of the error fnorm(est-tnsr)
all_resids
vector containing the Frobenius norm of error for all the iterations
The length of ranks
must match tnsr@num_modes-1
.
H. Lu, K. Plataniotis, A. Venetsanopoulos, "Mpca: Multilinear principal component analysis of tensor objects". IEEE Trans. Neural networks, 2008.
### How to retrieve faces_tnsr from figshare # faces_tnsr <- load_orl() # subject <- faces_tnsr[,,21,] dummy_faces_tnsr <- rand_tensor(c(92,112,40,10)) subject <- dummy_faces_tnsr[,,21,] mpcaD <- mpca(subject, ranks=c(10, 10)) mpcaD$conv mpcaD$norm_percent plot(mpcaD$all_resids)
### How to retrieve faces_tnsr from figshare # faces_tnsr <- load_orl() # subject <- faces_tnsr[,,21,] dummy_faces_tnsr <- rand_tensor(c(92,112,40,10)) subject <- dummy_faces_tnsr[,,21,] mpcaD <- mpca(subject, ranks=c(10, 10)) mpcaD$conv mpcaD$norm_percent plot(mpcaD$all_resids)
Overloads elementwise operators for tensors, arrays, and vectors that are conformable (have the same modes).
## S4 method for signature 'Tensor,Tensor' Ops(e1, e2)
## S4 method for signature 'Tensor,Tensor' Ops(e1, e2)
e1 |
left-hand object |
e2 |
right-hand object |
tnsr <- rand_tensor(c(3,4,5)) tnsr2 <- rand_tensor(c(3,4,5)) tnsrsum <- tnsr + tnsr2 tnsrdiff <- tnsr - tnsr2 tnsrelemprod <- tnsr * tnsr2 tnsrelemquot <- tnsr / tnsr2 for (i in 1:3L){ for (j in 1:4L){ for (k in 1:5L){ stopifnot(tnsrsum@data[i,j,k]==tnsr@data[i,j,k]+tnsr2@data[i,j,k]) stopifnot(tnsrdiff@data[i,j,k]==(tnsr@data[i,j,k]-tnsr2@data[i,j,k])) stopifnot(tnsrelemprod@data[i,j,k]==tnsr@data[i,j,k]*tnsr2@data[i,j,k]) stopifnot(tnsrelemquot@data[i,j,k]==tnsr@data[i,j,k]/tnsr2@data[i,j,k]) } } }
tnsr <- rand_tensor(c(3,4,5)) tnsr2 <- rand_tensor(c(3,4,5)) tnsrsum <- tnsr + tnsr2 tnsrdiff <- tnsr - tnsr2 tnsrelemprod <- tnsr * tnsr2 tnsrelemquot <- tnsr / tnsr2 for (i in 1:3L){ for (j in 1:4L){ for (k in 1:5L){ stopifnot(tnsrsum@data[i,j,k]==tnsr@data[i,j,k]+tnsr2@data[i,j,k]) stopifnot(tnsrdiff@data[i,j,k]==(tnsr@data[i,j,k]-tnsr2@data[i,j,k])) stopifnot(tnsrelemprod@data[i,j,k]==tnsr@data[i,j,k]*tnsr2@data[i,j,k]) stopifnot(tnsrelemquot@data[i,j,k]==tnsr@data[i,j,k]/tnsr2@data[i,j,k]) } } }
A wrapper function to image() to allow easy visualization of faces_tnsr, the ORL Face Dataset. The data is now stored in figshare https://ndownloader.figshare.com/files/28005669
plot_orl(subject = 1, condition = 1)
plot_orl(subject = 1, condition = 1)
subject |
which subject to plot (1-40) |
condition |
which lighting condition (1-10) |
AT&T Laboratories Cambridge. https://www.kaggle.com/kasikrit/att-database-of-faces
F. Samaria, A. Harter, "Parameterisation of a Stochastic Model for Human Face Identification". IEEE Workshop on Applications of Computer Vision 1994.
Extend print for Tensor
## S4 method for signature 'Tensor' print(x, ...)
## S4 method for signature 'Tensor' print(x, ...)
x |
the Tensor instance |
... |
additional parameters to be passed into print() |
print(x,...)
tnsr <- rand_tensor() print(tnsr)
tnsr <- rand_tensor() print(tnsr)
The default Population Value Decomposition (PVD) of a series of 2D images. Constructs population-level matrices P, V, and D to account for variances within as well as across the images. Structurally similar to Tucker (tucker
) and GLRAM (mpca
), but retains crucial differences. Requires 2*n3 + 2
parameters to specified the final ranks of P, V, and D, where n3 is the third mode (how many images are in the set). Consult Crainiceanu et al. (2013) for the construction and rationale behind the PVD model.
pvd(tnsr, uranks = NULL, wranks = NULL, a = NULL, b = NULL)
pvd(tnsr, uranks = NULL, wranks = NULL, a = NULL, b = NULL)
tnsr |
3-Tensor with the third mode being the measurement mode |
uranks |
ranks of the U matrices |
wranks |
ranks of the W matrices |
a |
rank of |
b |
rank of |
The PVD is not an iterative method, but instead relies on n3 + 2
separate PCA decompositions. The third mode is for how many images are in the set.
a list containing the following:
P
population-level matrix P = U%*%t(U)
, where U is constructed by stacking the truncated left eigenvectors of slicewise PCA along the third mode
V
a list of image-level core matrices
D
population-leve matrix D = W%*%t(W)
, where W is constructed by stacking the truncated right eigenvectors of slicewise PCA along the third mode
est
estimate of tnsr
after compression
norm_percent
the percent of Frobenius norm explained by the approximation
fnorm_resid
the Frobenius norm of the error fnorm(est-tnsr)
C. Crainiceanu, B. Caffo, S. Luo, V. Zipunnikov, N. Punjabi, "Population value decomposition: a framework for the analysis of image populations". Journal of the American Statistical Association, 2013.
### How to retrieve faces_tnsr from figshare # faces_tnsr <- load_orl() # subject <- faces_tnsr[,,8,] dummy_faces_tnsr <- rand_tensor(c(92,112,40,10)) subject <- dummy_faces_tnsr[,,8,] pvdD <- pvd(subject, uranks=rep(46,10), wranks=rep(56,10), a=46, b=56) plot(pvdD$fnorm_resid)
### How to retrieve faces_tnsr from figshare # faces_tnsr <- load_orl() # subject <- faces_tnsr[,,8,] dummy_faces_tnsr <- rand_tensor(c(92,112,40,10)) subject <- dummy_faces_tnsr[,,8,] pvdD <- pvd(subject, uranks=rep(46,10), wranks=rep(56,10), a=46, b=56) plot(pvdD$fnorm_resid)
Generate a Tensor with specified modes with iid normal(0,1) entries.
rand_tensor(modes = c(3, 4, 5), drop = FALSE)
rand_tensor(modes = c(3, 4, 5), drop = FALSE)
modes |
the modes of the output Tensor |
drop |
whether or not modes equal to 1 should be dropped |
a Tensor object with modes given by modes
Default rand_tensor()
generates a 3-Tensor with modes c(3,4,5)
.
rand_tensor() rand_tensor(c(4,4,4)) rand_tensor(c(10,2,1),TRUE)
rand_tensor() rand_tensor(c(4,4,4)) rand_tensor(c(10,2,1),TRUE)
DEPRECATED. Please see k_fold
.
rs_fold(mat, m = NULL, modes = NULL)
rs_fold(mat, m = NULL, modes = NULL)
mat |
matrix to be folded |
m |
the mode corresponding to rs_unfold |
modes |
the original modes of the tensor |
DEPRECATED. Please see k_unfold-methods
and unfold-methods
.
rs_unfold(tnsr, m) ## S4 method for signature 'Tensor' rs_unfold(tnsr, m = NULL)
rs_unfold(tnsr, m) ## S4 method for signature 'Tensor' rs_unfold(tnsr, m = NULL)
tnsr |
Tensor instance |
m |
mode to be unfolded on |
rs_unfold(tnsr,m=NULL)
Extend show for Tensor
## S4 method for signature 'Tensor' show(object)
## S4 method for signature 'Tensor' show(object)
object |
the Tensor instance |
show(object)
tnsr <- rand_tensor() tnsr
tnsr <- rand_tensor() tnsr
Implements T-MULT based on block circulant matrices (Kilmer et al. 2013) for 3-tensors.
t_mult(x, y)
t_mult(x, y)
x |
a 3-tensor |
y |
another 3-tensor |
Uses the Fast Fourier Transform (FFT) speed up suggested by Kilmer et al. 2013 instead of explicitly constructing the block circulant matrix. For the mathematical details of T-MULT, see Kilmer et al. (2013).
tensor product between x
and y
This only works (so far) between 3-Tensors.
M. Kilmer, K. Braman, N. Hao, and R. Hoover, "Third-order tensors as operators on matrices: a theoretical and computational framework with applications in imaging". SIAM Journal on Matrix Analysis and Applications 2013.
tnsr <- new("Tensor",3L,c(3L,4L,5L),data=runif(60)) tnsr2 <- new("Tensor",3L,c(4L,3L,5L),data=runif(60)) t_mult(tnsr, tnsr2)
tnsr <- new("Tensor",3L,c(3L,4L,5L),data=runif(60)) tnsr2 <- new("Tensor",3L,c(4L,3L,5L),data=runif(60)) t_mult(tnsr, tnsr2)
TSVD for a 3-Tensor. Constructs 3-Tensors U, S, V
such that tnsr = t_mult(t_mult(U,S),t(V))
. U
and V
are orthgonal 3-Tensors with orthogonality defined in Kilmer et al. (2013), and S
is a 3-Tensor consists of facewise diagonal matrices. For more details on the TSVD, consult Kilmer et al. (2013).
t_svd(tnsr)
t_svd(tnsr)
tnsr |
3-Tensor to decompose via TSVD |
a list containing the following:
U
the left orthgonal 3-Tensor
V
the right orthgonal 3-Tensor
S
the middle 3-Tensor consisting of face-wise diagonal matrices
Computation involves complex values, but if the inputs are real, then the outputs are also real. Some loss of precision occurs in the truncation of the imaginary components during the FFT and inverse FFT.
M. Kilmer, K. Braman, N. Hao, and R. Hoover, "Third-order tensors as operators on matrices: a theoretical and computational framework with applications in imaging". SIAM Journal on Matrix Analysis and Applications 2013.
tnsr <- rand_tensor() tsvdD <- t_svd(tnsr)
tnsr <- rand_tensor() tsvdD <- t_svd(tnsr)
Reconstruct the original 3-Tensor after it has been decomposed into U, S, V
via t_svd
.
t_svd_reconstruct(L)
t_svd_reconstruct(L)
L |
list that is an output from |
a 3-Tensor
tnsr <- rand_tensor(c(10,10,10)) tsvdD <- t_svd(tnsr) 1 - fnorm(t_svd_reconstruct(tsvdD)-tnsr)/fnorm(tnsr)
tnsr <- rand_tensor(c(10,10,10)) tsvdD <- t_svd(tnsr) 1 - fnorm(t_svd_reconstruct(tsvdD)-tnsr)/fnorm(tnsr)
Implements the tensor transpose based on block circulant matrices (Kilmer et al. 2013) for 3-tensors.
## S4 method for signature 'Tensor' t(x)
## S4 method for signature 'Tensor' t(x)
x |
a 3-tensor |
t(x)
tensor transpose of x
M. Kilmer, K. Braman, N. Hao, and R. Hoover, "Third-order tensors as operators on matrices: a theoretical and computational framework with applications in imaging". SIAM Journal on Matrix Analysis and Applications 2013.
tnsr <- rand_tensor() identical(t(tnsr)@data[,,1],t(tnsr@data[,,1])) identical(t(tnsr)@data[,,2],t(tnsr@data[,,5])) identical(t(t(tnsr)),tnsr)
tnsr <- rand_tensor() identical(t(tnsr)@data[,,1],t(tnsr@data[,,1])) identical(t(tnsr)@data[,,2],t(tnsr@data[,,5])) identical(t(t(tnsr)),tnsr)
Extend tail for Tensor
## S4 method for signature 'Tensor' tail(x, ...)
## S4 method for signature 'Tensor' tail(x, ...)
x |
the Tensor instance |
... |
additional parameters to be passed into tail() |
tail(x,...)
tnsr <- rand_tensor() tail(tnsr)
tnsr <- rand_tensor() tail(tnsr)
An S4 class for a tensor with arbitrary number of modes. The Tensor class extends the base 'array' class to include additional tensor manipulation (folding, unfolding, reshaping, subsetting) as well as a formal class definition that enables more explicit tensor algebra.
This can be seen as a wrapper class to the base array
class. While it is possible to create an instance using new
, it is also possible to do so by passing the data into as.tensor
.
Each slot of a Tensor instance can be obtained using @
.
The following methods are overloaded for the Tensor class: dim-methods
, head-methods
, tail-methods
, print-methods
, show-methods
, element-wise array operations, array subsetting (extract via ‘[’), array subset replacing (replace via ‘[<-’), and tperm-methods
, which is a wrapper around the base aperm
method.
To sum across any one mode of a tenor, use the function modeSum-methods
. To compute the mean across any one mode, use modeMean-methods
.
You can always unfold any Tensor into a matrix, and the unfold-methods
, k_unfold-methods
, and matvec-methods
methods are for that purpose. The output can be kept as a Tensor with 2 modes or a matrix
object. The vectorization function is also provided as vec
. See the attached vignette for a visualization of the different unfoldings.
Conversion from array
/matrix
to Tensor is facilitated via as.tensor
. To convert from a Tensor instance, simply invoke @data
.
The Frobenius norm of the Tensor is given by fnorm-methods
, while the inner product between two Tensors (of equal modes) is given by innerProd-methods
. You can also sum through any one mode to obtain the K-1 Tensor sum using modeSum-methods
. modeMean-methods
provides similar functionality to obtain the K-1 Tensor mean. These are primarily meant to be used internally but may be useful in doing statistics with Tensors.
For Tensors with 3 modes, we also overloaded t
(transpose) defined by Kilmer et.al (2013). See t-methods
.
To create a Tensor with i.i.d. random normal(0, 1) entries, see rand_tensor
.
number of modes (integer)
vector of modes (integer), aka sizes/extents/dimensions
actual data of the tensor, which can be 'array' or 'vector'
signature(tnsr = "Tensor")
: ...
signature(tnsr = "Tensor")
: ...
signature(tnsr = "Tensor")
: ...
signature(tnsr = "Tensor")
: ...
signature(tnsr = "Tensor")
: ...
signature(tnsr = "Tensor")
: ...
signature(.Object = "Tensor")
: ...
signature(tnsr1 = "Tensor", tnsr2 = "Tensor")
: ...
signature(tnsr = "Tensor")
: ...
signature(tnsr = "Tensor")
: ...
signature(e1 = "array", e2 = "Tensor")
: ...
signature(e1 = "numeric", e2 = "Tensor")
: ...
signature(e1 = "Tensor", e2 = "array")
: ...
signature(e1 = "Tensor", e2 = "numeric")
: ...
signature(e1 = "Tensor", e2 = "Tensor")
: ...
signature(tnsr = "Tensor")
: ...
signature(tnsr = "Tensor")
: ...
signature(tnsr = "Tensor")
: ...
signature(tnsr = "Tensor")
: ...
signature(tnsr = "Tensor")
: ...
signature(tnsr = "Tensor")
: ...
signature(tnsr = "Tensor")
: ...
signature(tnsr = "Tensor")
: ...
All of the decompositions and regression models in this package require a Tensor input.
James Li [email protected]
James Li, Jacob Bien, Martin T. Wells (2018). rTensor: An R Package for Multidimensional Array (Tensor) Unfolding, Multiplication, and Decomposition. Journal of Statistical Software, 87(10), 1-31. URL http://www.jstatsoft.org/v087/i10/.
tnsr <- rand_tensor() class(tnsr) tnsr print(tnsr) dim(tnsr) tnsr@num_modes tnsr@data
tnsr <- rand_tensor() class(tnsr) tnsr print(tnsr) dim(tnsr) tnsr@num_modes tnsr@data
Overloads aperm
for Tensor class for convenience.
tperm(tnsr, perm, ...) ## S4 method for signature 'Tensor' tperm(tnsr, perm, ...)
tperm(tnsr, perm, ...) ## S4 method for signature 'Tensor' tperm(tnsr, perm, ...)
tnsr |
the Tensor instance |
perm |
the new permutation of the current modes |
... |
additional parameters to be passed into |
tperm(tnsr,perm=NULL,...)
tnsr <- rand_tensor(c(3,4,5)) dim(tperm(tnsr,perm=c(2,1,3))) dim(tperm(tnsr,perm=c(1,3,2)))
tnsr <- rand_tensor(c(3,4,5)) dim(tperm(tnsr,perm=c(2,1,3))) dim(tperm(tnsr,perm=c(1,3,2)))
Contracted (m-Mode) product between a Tensor of arbitrary number of modes and a list of matrices. The result is folded back into Tensor.
ttl(tnsr, list_mat, ms = NULL)
ttl(tnsr, list_mat, ms = NULL)
tnsr |
Tensor object with K modes |
list_mat |
a list of matrices |
ms |
a vector of modes to contract on (order should match the order of |
Performs ttm
repeated for a single Tensor and a list of matrices on multiple modes. For instance, suppose we want to do multiply a Tensor object tnsr
with three matrices mat1
, mat2
, mat3
on modes 1, 2, and 3. We could do ttm(ttm(ttm(tnsr,mat1,1),mat2,2),3)
, or we could do ttl(tnsr,list(mat1,mat2,mat3),c(1,2,3))
. The order of the matrices in the list should obviously match the order of the modes. This is a common operation for various Tensor decompositions such as CP and Tucker. For the math on the m-Mode Product, see Kolda and Bader (2009).
Tensor object with K modes
The returned Tensor does not drop any modes equal to 1.
T. Kolda, B. Bader, "Tensor decomposition and applications". SIAM Applied Mathematics and Applications 2009.
tnsr <- new("Tensor",3L,c(3L,4L,5L),data=runif(60)) lizt <- list('mat1' = matrix(runif(30),ncol=3), 'mat2' = matrix(runif(40),ncol=4), 'mat3' = matrix(runif(50),ncol=5)) ttl(tnsr,lizt,ms=c(1,2,3))
tnsr <- new("Tensor",3L,c(3L,4L,5L),data=runif(60)) lizt <- list('mat1' = matrix(runif(30),ncol=3), 'mat2' = matrix(runif(40),ncol=4), 'mat3' = matrix(runif(50),ncol=5)) ttl(tnsr,lizt,ms=c(1,2,3))
Contracted (m-Mode) product between a Tensor of arbitrary number of modes and a matrix. The result is folded back into Tensor.
ttm(tnsr, mat, m = NULL)
ttm(tnsr, mat, m = NULL)
tnsr |
Tensor object with K modes |
mat |
input matrix with same number columns as the |
m |
the mode to contract on |
By definition, rs_unfold(ttm(tnsr,mat),m) = mat%*%rs_unfold(tnsr,m)
, so the number of columns in mat
must match the m
th mode of tnsr
. For the math on the m-Mode Product, see Kolda and Bader (2009).
a Tensor object with K modes
The m
th mode of tnsr
must match the number of columns in mat
. By default, the returned Tensor does not drop any modes equal to 1.
T. Kolda, B. Bader, "Tensor decomposition and applications". SIAM Applied Mathematics and Applications 2009.
tnsr <- new("Tensor",3L,c(3L,4L,5L),data=runif(60)) mat <- matrix(runif(50),ncol=5) ttm(tnsr,mat,m=3)
tnsr <- new("Tensor",3L,c(3L,4L,5L),data=runif(60)) mat <- matrix(runif(50),ncol=5) ttm(tnsr,mat,m=3)
The Tucker decomposition of a tensor. Approximates a K-Tensor using a n-mode product of a core tensor (with modes specified by ranks
) with orthogonal factor matrices. If there is no truncation in one of the modes, then this is the same as the MPCA, mpca
. If there is no truncation in all the modes (i.e. ranks = tnsr@modes
), then this is the same as the HOSVD, hosvd
. This is an iterative algorithm, with two possible stopping conditions: either relative error in Frobenius norm has gotten below tol
, or the max_iter
number of iterations has been reached. For more details on the Tucker decomposition, consult Kolda and Bader (2009).
tucker(tnsr, ranks = NULL, max_iter = 25, tol = 1e-05)
tucker(tnsr, ranks = NULL, max_iter = 25, tol = 1e-05)
tnsr |
Tensor with K modes |
ranks |
a vector of the modes of the output core Tensor |
max_iter |
maximum number of iterations if error stays above |
tol |
relative Frobenius norm error tolerance |
Uses the Alternating Least Squares (ALS) estimation procedure also known as Higher-Order Orthogonal Iteration (HOOI). Intialized using a (Truncated-)HOSVD. A progress bar is included to help monitor operations on large tensors.
a list containing the following:
Z
the core tensor, with modes specified by ranks
U
a list of orthgonal factor matrices - one for each mode, with the number of columns of the matrices given by ranks
conv
whether or not resid
< tol
by the last iteration
est
estimate of tnsr
after compression
norm_percent
the percent of Frobenius norm explained by the approximation
fnorm_resid
the Frobenius norm of the error fnorm(est-tnsr)
all_resids
vector containing the Frobenius norm of error for all the iterations
The length of ranks
must match tnsr@num_modes
.
T. Kolda, B. Bader, "Tensor decomposition and applications". SIAM Applied Mathematics and Applications 2009.
tnsr <- rand_tensor(c(4,4,4,4)) tuckerD <- tucker(tnsr,ranks=c(2,2,2,2)) tuckerD$conv tuckerD$norm_percent plot(tuckerD$all_resids)
tnsr <- rand_tensor(c(4,4,4,4)) tuckerD <- tucker(tnsr,ranks=c(2,2,2,2)) tuckerD$conv tuckerD$norm_percent plot(tuckerD$all_resids)
Unfolds the tensor into a matrix, with the modes in rs
onto the rows and modes in cs
onto the columns. Note that c(rs,cs)
must have the same elements (order doesn't matter) as x@modes
. Within the rows and columns, the order of the unfolding is determined by the order of the modes. This convention is consistent with Kolda and Bader (2009).
unfold(tnsr, row_idx, col_idx) ## S4 method for signature 'Tensor' unfold(tnsr, row_idx = NULL, col_idx = NULL)
unfold(tnsr, row_idx, col_idx) ## S4 method for signature 'Tensor' unfold(tnsr, row_idx = NULL, col_idx = NULL)
tnsr |
the Tensor instance |
row_idx |
the indices of the modes to map onto the row space |
col_idx |
the indices of the modes to map onto the column space |
For Row Space Unfolding or m-mode Unfolding, see rs_unfold-methods
. For Column Space Unfolding or matvec, see cs_unfold-methods
.
vec-methods
returns the vectorization of the tensor.
unfold(tnsr,row_idx=NULL,col_idx=NULL)
matrix with prod(row_idx)
rows and prod(col_idx)
columns
T. Kolda, B. Bader, "Tensor decomposition and applications". SIAM Applied Mathematics and Applications 2009.
k_unfold-methods
and matvec-methods
tnsr <- rand_tensor() matT3<-unfold(tnsr,row_idx=2,col_idx=c(3,1))
tnsr <- rand_tensor() matT3<-unfold(tnsr,row_idx=2,col_idx=c(3,1))
The inverse operation to matvec-methods
, turning a matrix into a Tensor. For a full account of matrix folding/unfolding operations, consult Kolda and Bader (2009).
unmatvec(mat, modes = NULL)
unmatvec(mat, modes = NULL)
mat |
matrix to be folded into a Tensor |
modes |
the modes of the output Tensor |
Tensor object with modes given by modes
T. Kolda, B. Bader, "Tensor decomposition and applications". SIAM Applied Mathematics and Applications 2009.
tnsr <- new("Tensor",3L,c(3L,4L,5L),data=runif(60)) matT1<-matvec(tnsr) identical(unmatvec(matT1,modes=c(3,4,5)),tnsr)
tnsr <- new("Tensor",3L,c(3L,4L,5L),data=runif(60)) matT1<-matvec(tnsr) identical(unmatvec(matT1,modes=c(3,4,5)),tnsr)
Turns the tensor into a single vector, following the convention that earlier indices vary slower than later indices.
vec(tnsr) ## S4 method for signature 'Tensor' vec(tnsr)
vec(tnsr) ## S4 method for signature 'Tensor' vec(tnsr)
tnsr |
the Tensor instance |
vec(tnsr)
vector with length prod(x@modes)
T. Kolda, B. Bader, "Tensor decomposition and applications". SIAM Applied Mathematics and Applications 2009.
tnsr <- rand_tensor(c(4,5,6,7)) vec(tnsr)
tnsr <- rand_tensor(c(4,5,6,7)) vec(tnsr)