Title: | R package for sparse and out-of-core arithmetic and decomposition of Tensor |
---|---|
Description: | DelayedTensor operates Tensor arithmetic directly on DelayedArray object. DelayedTensor provides some generic function related to Tensor arithmetic/decompotision and dispatches it on the DelayedArray class. DelayedTensor also suppors Tensor contraction by einsum function, which is inspired by numpy einsum. |
Authors: | Koki Tsuyuzaki [aut, cre] |
Maintainer: | Koki Tsuyuzaki <[email protected]> |
License: | Artistic-2.0 |
Version: | 1.13.0 |
Built: | 2024-10-30 05:28:23 UTC |
Source: | https://github.com/bioc/DelayedTensor |
DelayedTensor operates Tensor arithmetic directly on DelayedArray object. DelayedTensor provides some generic function related to Tensor arithmetic/decompotision and dispatches it on the DelayedArray class. DelayedTensor also suppors Tensor contraction by einsum function, which is inspired by numpy einsum.
The DESCRIPTION file:
Package: | DelayedTensor |
Type: | Package |
Title: | R package for sparse and out-of-core arithmetic and decomposition of Tensor |
Version: | 1.13.0 |
Authors@R: | c(person("Koki", "Tsuyuzaki", role = c("aut", "cre"), email = "[email protected]")) |
Depends: | R (>= 4.1.0) |
Imports: | methods, utils, S4Arrays, SparseArray, DelayedArray (>= 0.31.8), HDF5Array, BiocSingular, rTensor, DelayedRandomArray (>= 1.13.1), irlba, Matrix, einsum, |
Suggests: | markdown, rmarkdown, BiocStyle, knitr, testthat, magrittr, dplyr, reticulate |
Description: | DelayedTensor operates Tensor arithmetic directly on DelayedArray object. DelayedTensor provides some generic function related to Tensor arithmetic/decompotision and dispatches it on the DelayedArray class. DelayedTensor also suppors Tensor contraction by einsum function, which is inspired by numpy einsum. |
License: | Artistic-2.0 |
biocViews: | Software, Infrastructure, DataRepresentation, DimensionReduction |
LazyData: | true |
LazyDataCompression: | xz |
VignetteBuilder: | knitr |
BugReports: | https://github.com/rikenbit/DelayedTensor/issues |
Repository: | https://bioc.r-universe.dev |
RemoteUrl: | https://github.com/bioc/DelayedTensor |
RemoteRef: | HEAD |
RemoteSha: | 35f745ff4ac57c3f7acba76b2d2deecebdd13c5f |
Author: | Koki Tsuyuzaki [aut, cre] |
Maintainer: | Koki Tsuyuzaki <[email protected]> |
Index of help topics:
DelayedDiagonalArray Diagonal DelayedArray DelayedTensor-package R package for sparse and out-of-core arithmetic and decomposition of Tensor cbind_list Mode-binding against list cp-methods Canonical Polyadic Decomposition cs_fold-methods Column Space Folding of 2D DelayedArray cs_unfold-methods Tensor Column Space Unfolding of DelayedArray diag-methods DelayedArray Diagonals einsum Einstein Summation of DelayedArray fnorm-methods Tensor Frobenius Norm of DelayedArray fold-methods Tensor folding of 2D DelayedArray getSparse Getter of the intermediate/output DelayedArray object in DelayedTensor getVerbose Getter function to control the verbose messages from DelayedTensor hadamard-methods Hadamard Product of DelayedArray hadamard_list Hadamard Product against list hosvd-methods (Truncated-)Higher-order SVD human_mid_brain Matrix object of human mid brain data innerProd-methods Tensors Inner Product of DelayedArray k_fold-methods k-mode Folding of 2D DelayedArray k_unfold-methods Tensor k-mode Unfolding of DelayedArray khatri_rao-methods Khatri-Rao Product of DelayedArray khatri_rao_list Khatri-Rao Product against list kronecker-methods Kronecker Product of DelayedArray kronecker_list Kronecker Product against list list_rep Replicate of arbitrary object matvec-methods Tensor Matvec Unfolding of DelayedArray modeMean-methods Tensor Mean Across Single Mode of DelayedArray modeSum-methods Tensor Sum Across Single Mode of DelayedArray modebind_list Mode-binding against list mouse_mid_brain Matrix object of mouse mid brain data mpca-methods Multilinear Principal Components Analysis outerProd-methods Tensors Outer Product of DelayedArray pvd-methods Population Value Decomposition rbind_list Mode-binding against list rs_fold-methods Row Space Folding of 2D DelayedArray rs_unfold-methods Tensor Row Space Unfolding of DelayedArray setSparse Setter to set the intermediate DelayedArray object in DelayedTensor setVerbose Setter to set the verbose mode of DelayedTensor ttl DelayedArray Times List ttm-methods Tensor Times Matrix (m-Mode Product) tucker-methods Tucker Decomposition unfold-methods Tensor Unfolding of 2D DelayedArray unmatvec-methods Unmatvec Folding of 2D DelayedArray vec-methods Tensor Vectorization of DelayedArray
Koki Tsuyuzaki [aut, cre]
Maintainer: Koki Tsuyuzaki <[email protected]>
# Unfold operations
unfold
, k_unfold
, matvec
,
rs_unfold
, cs_unfold
, ttl
# Fold operations
fold
, k_fold
, unmatvec
,
rs_fold
, cs_fold
, ttm
# Vectorization
vec
# Norm operations
fnorm
, innerProd
# Diagonal operations / Diagonal Tensor
diag
, DelayedDiagonalArray
# Mode-wise operations
modeSum
, modeMean
# Tensor product operations
hadamard
, hadamard_list
,
kronecker
, kronecker_list
,
khatri_rao
, khatri_rao_list
# Utilities
list_rep
, modebind_list
,
rbind_list
, cbind_list
# Decomposition operations
hosvd
, cp
, tucker
,
mpca
, pvd
# Einsum operation
einsum
ls("package:DelayedTensor")
ls("package:DelayedTensor")
Returns the binded DelayedArray in column space.
cbind_list(L)
cbind_list(L)
L |
list of 2D DelayedArray |
This is a wrapper function to modebind_list
,
when the DelayedArrays are 2D.
2D DelayedArray object
The dimensions of column in each DelayedArray must match.
library("DelayedRandomArray") dlizt <- list( 'darr1' = RandomUnifArray(c(2,3)), 'darr2' = RandomUnifArray(c(2,3))) cbind_list(dlizt)
library("DelayedRandomArray") dlizt <- list( 'darr1' = RandomUnifArray(c(2,3)), 'darr2' = RandomUnifArray(c(2,3))) cbind_list(dlizt)
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 *darr@num_modes
vectors in the output,
stored in darr@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(darr, num_components=NULL, max_iter=25, tol=1e-05) ## S4 method for signature 'DelayedArray' cp(darr, num_components, max_iter, tol)
cp(darr, num_components=NULL, max_iter=25, tol=1e-05) ## S4 method for signature 'DelayedArray' cp(darr, num_components, max_iter, tol)
darr |
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 |
This function is an extension of the cp
by DelayedArray.
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 darr
after compression
fnorm_resid
the Frobenius norm of the error
fnorm(est-darr)
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.
library("DelayedRandomArray") darr <- RandomUnifArray(c(3,4,5)) cp(darr, num_components=2)
library("DelayedRandomArray") darr <- RandomUnifArray(c(3,4,5)) cp(darr, num_components=2)
The inverse operation to cs_unfold
.
cs_fold(mat, m = NULL, modes = NULL) ## S4 method for signature 'DelayedArray' cs_fold(mat, m, modes)
cs_fold(mat, m = NULL, modes = NULL) ## S4 method for signature 'DelayedArray' cs_fold(mat, m, modes)
mat |
DelayedArray object (only 2D) |
m |
the mode corresponding to cs_unfold |
modes |
the original modes of the DelayedArray |
This function is an extension of the cs_fold
by DelayedArray.
This is a wrapper function to fold
.
DelayedArray (higher than 2D)
T. Kolda, B. Bader, "Tensor decomposition and applications". SIAM Applied Mathematics and Applications 2009.
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) matT3 <- DelayedTensor::cs_unfold(darr, m=3) identical( as.array(DelayedTensor::cs_fold(matT3, m=3, modes=c(2,3,4))), as.array(darr))
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) matT3 <- DelayedTensor::cs_unfold(darr, m=3) identical( as.array(DelayedTensor::cs_fold(matT3, m=3, modes=c(2,3,4))), as.array(darr))
cs_unfold(darr, m) ## S4 method for signature 'DelayedArray' cs_unfold(darr, m)
cs_unfold(darr, m) ## S4 method for signature 'DelayedArray' cs_unfold(darr, m)
darr |
DelayedArray object |
m |
mode to be unfolded on |
This function is an extension of the cs_unfold
by DelayedArray.
This is a wrapper function to unfold
.
DelayedArray (2D)
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) DelayedTensor::cs_unfold(darr, m=3)
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) DelayedTensor::cs_unfold(darr, m=3)
Constructor of the diagonal of a DelayedArray.
DelayedDiagonalArray(shape, value)
DelayedDiagonalArray(shape, value)
shape |
Shape of DelayedArray (mode of Tensor) |
value |
either a single value or a vector. This argument is optional. If nothing is specified, 1s are filled with each diagonal element. |
DelayedArray object
T. Kolda, B. Bader, "Tensor decomposition and applications". SIAM Applied Mathematics and Applications 2009.
darr <- DelayedDiagonalArray(2:4, 5) DelayedTensor::diag(darr)
darr <- DelayedDiagonalArray(2:4, 5) DelayedTensor::diag(darr)
Extract or replace the diagonal of a DelayedArray, or substitute the elements to the diagonal DelayedArray.
diag(darr) diag(darr) <- value ## S4 method for signature 'DelayedArray' diag(darr) ## S4 replacement method for signature 'DelayedArray' diag(darr) <- value
diag(darr) diag(darr) <- value ## S4 method for signature 'DelayedArray' diag(darr) ## S4 replacement method for signature 'DelayedArray' diag(darr) <- value
darr |
DelayedArray object |
value |
either a single value or a vector of length equal to that
of the current diagonal. Should be of a mode which can be coerced
to that of |
See also DelayedDiagonalArray
or diag
.
1D DelayedArray (vector) with length min(dim(darr))
T. Kolda, B. Bader, "Tensor decomposition and applications". SIAM Applied Mathematics and Applications 2009.
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) DelayedTensor::diag(darr) DelayedTensor::diag(darr)[1] <- 11111 DelayedTensor::diag(darr)[2] <- 22222 DelayedTensor::diag(darr)
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) DelayedTensor::diag(darr) DelayedTensor::diag(darr)[1] <- 11111 DelayedTensor::diag(darr)[2] <- 22222 DelayedTensor::diag(darr)
Einstein summation is a convenient and concise notation for operations on n-dimensional arrays.
NOTE: Sparse mode of einsum is not available for now.
einsum(subscripts, ...)
einsum(subscripts, ...)
subscripts |
a string in Einstein notation where arrays
are separated by ',' and the result is separated by '->'. For
example |
... |
the DelayedArrays that are combined. |
This function is an extension of the einsum
by DelayedArray.
The einsum
function returns an array with one dimension for each index
in the result of the subscripts
.
For example "ij,jk->ik"
produces a 2-dimensional array,
"abc,cd,de->abe"
produces a 3-dimensional array.
library("DelayedArray") library("DelayedRandomArray") darr1 <- RandomUnifArray(c(4,8)) darr2 <- RandomUnifArray(c(8,3)) # Matrix Multiply darr1 %*% darr2 DelayedTensor::einsum("ij,jk -> ik", darr1, darr2) # Diag mat_sq <- RandomUnifArray(c(4,4)) DelayedTensor::diag(mat_sq) einsum("ii->i", mat_sq) # Trace sum(DelayedTensor::diag(mat_sq)) einsum("ii->", mat_sq) # Scalar product darr3 <- RandomUnifArray(c(4,8)) darr3 * darr1 einsum("ij,ij->ij", darr3, darr1) # Transpose t(darr1) einsum("ij->ji", darr1) # Batched L2 norm arr1 <- as.array(darr1) arr3 <- as.array(darr3) darr4 <- DelayedArray(array(c(arr1, arr3), dim = c(dim(arr1), 2))) c(sum(darr1^2), sum(darr3^2)) einsum("ijb,ijb->b", darr4, darr4)
library("DelayedArray") library("DelayedRandomArray") darr1 <- RandomUnifArray(c(4,8)) darr2 <- RandomUnifArray(c(8,3)) # Matrix Multiply darr1 %*% darr2 DelayedTensor::einsum("ij,jk -> ik", darr1, darr2) # Diag mat_sq <- RandomUnifArray(c(4,4)) DelayedTensor::diag(mat_sq) einsum("ii->i", mat_sq) # Trace sum(DelayedTensor::diag(mat_sq)) einsum("ii->", mat_sq) # Scalar product darr3 <- RandomUnifArray(c(4,8)) darr3 * darr1 einsum("ij,ij->ij", darr3, darr1) # Transpose t(darr1) einsum("ij->ji", darr1) # Batched L2 norm arr1 <- as.array(darr1) arr3 <- as.array(darr3) darr4 <- DelayedArray(array(c(arr1, arr3), dim = c(dim(arr1), 2))) c(sum(darr1^2), sum(darr3^2)) einsum("ijb,ijb->b", darr4, darr4)
Returns the Frobenius norm of the Tensor instance.
fnorm(darr) ## S4 method for signature 'DelayedArray' fnorm(darr)
fnorm(darr) ## S4 method for signature 'DelayedArray' fnorm(darr)
darr |
DelayedArray object |
This function is an extension of the fnorm
by DelayedArray.
numeric Frobenius norm of darr
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) fnorm(darr)
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) fnorm(darr)
General folding of a 2D DelayedArray into a higher-order DelayedArray(Tensor).
This is designed to be the inverse function to
unfold
, 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) ## S4 method for signature 'DelayedArray' fold(mat, row_idx, col_idx, modes)
fold(mat, row_idx = NULL, col_idx = NULL, modes = NULL) ## S4 method for signature 'DelayedArray' fold(mat, row_idx, col_idx, modes)
mat |
DelayedArray object (only 2D) |
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 DelayedArray |
This function is an extension of the fold
by DelayedArray.
DelayedArray object with modes given by modes
T. Kolda, B. Bader, "Tensor decomposition and applications". SIAM Applied Mathematics and Applications 2009.
unfold
, k_fold
,
unmatvec
,
rs_fold
, cs_fold
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) matT3 <- DelayedTensor::unfold(darr, row_idx=2, col_idx=c(3,1)) identical( as.array(DelayedTensor::fold(matT3, row_idx=2,col_idx=c(3,1), modes=c(2,3,4))), as.array(darr))
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) matT3 <- DelayedTensor::unfold(darr, row_idx=2, col_idx=c(3,1)) identical( as.array(DelayedTensor::fold(matT3, row_idx=2,col_idx=c(3,1), modes=c(2,3,4))), as.array(darr))
Whether the intermediate and output DelayedArray used in DelayedTensor is used as sparse tensor or not.
NOTE: Sparse mode is experimental! Whether it contributes to higher speed and lower memory is quite dependent on the sparsity of the DelayedArray, and the current implementation does not recognize the block size, which may cause Out-of-Memory errors.
getSparse()
getSparse()
TRUE or FALSE (Default: FALSE)
getSparse()
getSparse()
Returns the verbose setting of DelayedTensor functions.
getVerbose()
getVerbose()
TRUE or FALSE (Default: FALSE)
getVerbose()
getVerbose()
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 DelayedArray |
This function is an extension of the hadamard_list
by DelayedArray.
matrix that is the Hadamard product
The modes/dimensions of each element in the list must match.
khatri_rao
, khatri_rao_list
,
kronecker
, kronecker_list
,
hadamard
library("DelayedRandomArray") dlizt <- list( 'darr1' = RandomUnifArray(c(2,3,4)), 'darr2' = RandomUnifArray(c(2,3,4))) hadamard_list(dlizt)
library("DelayedRandomArray") dlizt <- list( 'darr1' = RandomUnifArray(c(2,3,4)), 'darr2' = RandomUnifArray(c(2,3,4))) hadamard_list(dlizt)
Returns the Hadamard product of two Tensors. Commonly used for n-mode products and various Tensor decompositions.
hadamard(darr1, darr2) ## S4 method for signature 'DelayedArray,DelayedArray' hadamard(darr1, darr2)
hadamard(darr1, darr2) ## S4 method for signature 'DelayedArray,DelayedArray' hadamard(darr1, darr2)
darr1 |
first DelayedArray object |
darr2 |
second DelayedArray object |
matrix that is the Hadamard product
The modes/dimensions of each element of two Tensors must match.
khatri_rao
, khatri_rao_list
,
kronecker
, kronecker_list
,
hadamard_list
library("DelayedRandomArray") darr1 <- RandomUnifArray(c(2,4)) darr2 <- RandomUnifArray(c(2,4)) hadamard(darr1, darr1)
library("DelayedRandomArray") darr1 <- RandomUnifArray(c(2,4)) darr2 <- RandomUnifArray(c(2,4)) hadamard(darr1, darr1)
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(darr, ranks=NULL) ## S4 method for signature 'DelayedArray' hosvd(darr, ranks)
hosvd(darr, ranks=NULL) ## S4 method for signature 'DelayedArray' hosvd(darr, ranks)
darr |
Tensor with K modes |
ranks |
a vector of desired modes in the output core tensor,
default is |
This function is an extension of the hosvd
by DelayedArray.
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 darr
after compression
fnorm_resid
the Frobenius norm of the error
fnorm(est-darr)
- if there was no truncation,
then this is on the order of mach_eps * fnorm.
The length of ranks
must match darr@num_modes
.
L. Lathauwer, B.Moor, J. Vanderwalle "A multilinear singular value decomposition". Journal of Matrix Analysis and Applications 2000.
library("DelayedRandomArray") darr <- RandomUnifArray(c(3,4,5)) hosvd(darr, ranks=c(2,1,3))
library("DelayedRandomArray") darr <- RandomUnifArray(c(3,4,5)) hosvd(darr, ranks=c(2,1,3))
A matrix with 500 rows (genes) * 1977 columns (cells).
data(human_mid_brain)
data(human_mid_brain)
The data matrix is downloaded from GEO Series GSE76381 (https://www.ncbi.nlm.nih.gov/geo/download/?acc=GSE76381&format=file). For the details, see inst/script/make-data.R.
Y-h. Taguchi and T. Turki (2019) Tensor Decomposition-Based Unsupervised Feature Extraction Applied to Single-Cell Gene Expression Analysis. Frontiers in Genetics, 10(864): 10:3389/fgene.2019.00864
data(human_mid_brain)
data(human_mid_brain)
Returns the inner product between two Tensors
innerProd(darr1, darr2) ## S4 method for signature 'DelayedArray,DelayedArray' innerProd(darr1, darr2)
innerProd(darr1, darr2) ## S4 method for signature 'DelayedArray,DelayedArray' innerProd(darr1, darr2)
darr1 |
first DelayedArray object |
darr2 |
second DelayedArray object |
This function is an extension of the innerProd
by DelayedArray.
inner product between darr1
and darr2
library("DelayedRandomArray") darr1 <- RandomUnifArray(c(2,3,4)) darr2 <- RandomUnifArray(c(2,3,4)) innerProd(darr1, darr2)
library("DelayedRandomArray") darr1 <- RandomUnifArray(c(2,3,4)) darr2 <- RandomUnifArray(c(2,3,4)) innerProd(darr1, darr2)
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(darr, m), m, dim(darr))
will result in the original Tensor.
k_fold(mat, m = NULL, modes = NULL) ## S4 method for signature 'DelayedArray' k_fold(mat, m, modes)
k_fold(mat, m = NULL, modes = NULL) ## S4 method for signature 'DelayedArray' k_fold(mat, m, modes)
mat |
DelayedArray object (only 2D) |
m |
the index of the mode that is mapped onto the row indices |
modes |
the modes of the output DelayedArray |
This function is an extension of the k_fold
by DelayedArray.
This is a wrapper function to fold
.
DelayedArray object with modes given by modes
T. Kolda, B. Bader, "Tensor decomposition and applications". SIAM Applied Mathematics and Applications 2009.
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) matT2 <- k_unfold(darr, m=2) identical( as.array(k_fold(matT2, m=2, modes=c(2,3,4))), as.array(darr))
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) matT2 <- k_unfold(darr, m=2) identical( as.array(k_fold(matT2, m=2, modes=c(2,3,4))), as.array(darr))
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(darr, m) ## S4 method for signature 'DelayedArray' k_unfold(darr, m)
k_unfold(darr, m) ## S4 method for signature 'DelayedArray' k_unfold(darr, m)
darr |
DelayedArray object |
m |
the index of the mode to unfold on |
This function is an extension of the k_unfold
by DelayedArray.
This is a wrapper function to unfold
.
See also k_unfold(darr, m=NULL)
matrix with dim(darr)[m]
rows and prod(dim(darr)[-m])
columns
T. Kolda, B. Bader, "Tensor decomposition and applications". SIAM Applied Mathematics and Applications 2009.
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) rs_unfold(darr, m=2)
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) rs_unfold(darr, m=2)
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 DelayedArray |
reverse |
whether or not to reverse the order |
This function is an extension of the khatri_rao_list
by DelayedArray.
matrix that is the Khatri-Rao product
The number of columns must match in every element of the input list.
hadamard
, hadamard_list
,
kronecker
, kronecker_list
,
khatri_rao
library("DelayedRandomArray") darr1 <- RandomUnifArray(c(2,3)) dlizt <- list( 'darr1' = RandomUnifArray(c(2,4)), 'darr2' = RandomUnifArray(c(3,4))) khatri_rao_list(dlizt)
library("DelayedRandomArray") darr1 <- RandomUnifArray(c(2,3)) dlizt <- list( 'darr1' = RandomUnifArray(c(2,4)), 'darr2' = RandomUnifArray(c(3,4))) khatri_rao_list(dlizt)
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(darr1, darr2) ## S4 method for signature 'DelayedArray,DelayedArray' khatri_rao(darr1, darr2)
khatri_rao(darr1, darr2) ## S4 method for signature 'DelayedArray,DelayedArray' khatri_rao(darr1, darr2)
darr1 |
first DelayedArray object |
darr2 |
second DelayedArray object |
This function is an extension of the khatri_rao
by DelayedArray.
matrix that is the Khatri-Rao product
The number of columns must match in the two inputs.
hadamard
, hadamard_list
,
kronecker
, kronecker_list
,
khatri_rao_list
library("DelayedRandomArray") darr1 <- RandomUnifArray(c(2,4)) darr2 <- RandomUnifArray(c(3,4)) khatri_rao(darr1, darr2)
library("DelayedRandomArray") darr1 <- RandomUnifArray(c(2,4)) darr2 <- RandomUnifArray(c(3,4)) khatri_rao(darr1, darr2)
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 DelayedArray |
This function is an extension of the kronecker_list
by DelayedArray.
matrix that is the Kronecker product
khatri_rao
, khatri_rao_list
,
hadamard
, hadamard_list
,
kronecker
library("DelayedRandomArray") dlizt <- list( 'darr1' = RandomUnifArray(c(2,3,4)), 'darr2' = RandomUnifArray(c(2,3,4))) kronecker_list(dlizt)
library("DelayedRandomArray") dlizt <- list( 'darr1' = RandomUnifArray(c(2,3,4)), 'darr2' = RandomUnifArray(c(2,3,4))) kronecker_list(dlizt)
Returns the Kronecker product of two Tensors. Commonly used for n-mode products and various Tensor decompositions.
kronecker(darr1, darr2) ## S4 method for signature 'DelayedArray,DelayedArray' kronecker(darr1, darr2)
kronecker(darr1, darr2) ## S4 method for signature 'DelayedArray,DelayedArray' kronecker(darr1, darr2)
darr1 |
first DelayedArray object |
darr2 |
second DelayedArray object |
matrix that is the Kronecker product
khatri_rao
, khatri_rao_list
,
hadamard
, hadamard_list
,
kronecker_list
library("DelayedRandomArray") darr1 <- RandomUnifArray(c(2,3)) darr2 <- RandomUnifArray(c(4,5)) kronecker(darr1, darr2)
library("DelayedRandomArray") darr1 <- RandomUnifArray(c(2,3)) darr2 <- RandomUnifArray(c(4,5)) kronecker(darr1, darr2)
Returns the replicates of base obejct x
.
list_rep(x, n=NULL)
list_rep(x, n=NULL)
x |
Any object |
n |
Number of replicate |
List
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) list_rep(darr, 3)
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) list_rep(darr, 3)
For 3-tensors only. Stacks the slices along the third mode.
matvec(darr) ## S4 method for signature 'DelayedArray' matvec(darr)
matvec(darr) ## S4 method for signature 'DelayedArray' matvec(darr)
darr |
DelayedArray object |
This function is an extension of the matvec
by DelayedArray.
This is a wrapper function to unfold
.
matrix with prod(dim(darr)[-m])
rows and dim(darr)[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.
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) matvec(darr)
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) matvec(darr)
Returns the binded DelayedArray in mode-m
.
modebind_list(L, m=NULL)
modebind_list(L, m=NULL)
L |
list of DelayedArray |
m |
list of DelayedArray |
DelayedArray object
The dimensions of mode m
must match.
library("DelayedRandomArray") dlizt <- list( 'darr1' = RandomUnifArray(c(2,3,4)), 'darr2' = RandomUnifArray(c(2,3,4))) modebind_list(dlizt, m=1) modebind_list(dlizt, m=2) modebind_list(dlizt, m=3)
library("DelayedRandomArray") dlizt <- list( 'darr1' = RandomUnifArray(c(2,3,4)), 'darr2' = RandomUnifArray(c(2,3,4))) modebind_list(dlizt, m=1) modebind_list(dlizt, m=2) modebind_list(dlizt, m=3)
Given a mode for a K-tensor, this returns the K-1 tensor resulting from taking the mean across that particular mode.
modeMean(darr, m = NULL, drop = FALSE) ## S4 method for signature 'DelayedArray' modeMean(darr, m, drop)
modeMean(darr, m = NULL, drop = FALSE) ## S4 method for signature 'DelayedArray' modeMean(darr, m, drop)
darr |
DelayedArray object |
m |
the index of the mode to average across |
drop |
whether or not mode m should be dropped |
This function is an extension of the modeMean
by DelayedArray.
NOTE: Sparse mode of modeMean is not available for now.
modeMean(darr, m=NULL, drop=FALSE)
K-1 or K Tensor, where K = length(dim(darr))
library("DelayedRandomArray") darr <- RandomUnifArray(c(1,2,3)) modeMean(darr, 1, drop=FALSE) modeMean(darr, 1, drop=TRUE) modeMean(darr, 2) modeMean(darr, 3)
library("DelayedRandomArray") darr <- RandomUnifArray(c(1,2,3)) modeMean(darr, 1, drop=FALSE) modeMean(darr, 1, drop=TRUE) modeMean(darr, 2) modeMean(darr, 3)
Given a mode for a K-tensor, this returns the K-1 tensor resulting from summing across that particular mode.
modeSum(darr, m = NULL, drop = FALSE) ## S4 method for signature 'DelayedArray' modeSum(darr, m, drop)
modeSum(darr, m = NULL, drop = FALSE) ## S4 method for signature 'DelayedArray' modeSum(darr, m, drop)
darr |
DelayedArray object |
m |
the index of the mode to sum across |
drop |
whether or not mode m should be dropped |
This function is an extension of the modeSum
by DelayedArray.
NOTE: Sparse mode of modeSum is not available for now.
modeSum(darr, m=NULL, drop=FALSE)
K-1 or K tensor, where K = length(dim(darr))
library("DelayedRandomArray") darr <- RandomUnifArray(c(1,2,3)) modeSum(darr, 1, drop=FALSE) modeSum(darr, 1, drop=TRUE) modeSum(darr, 2) modeSum(darr, 3)
library("DelayedRandomArray") darr <- RandomUnifArray(c(1,2,3)) modeSum(darr, 1, drop=FALSE) modeSum(darr, 1, drop=TRUE) modeSum(darr, 2) modeSum(darr, 3)
A matrix with 500 rows (genes) * 1907 columns (cells).
data(mouse_mid_brain)
data(mouse_mid_brain)
The data matrix is downloaded from GEO Series GSE76381 (https://www.ncbi.nlm.nih.gov/geo/download/?acc=GSE76381&format=file). For the details, see inst/script/make-data.R.
Y-h. Taguchi and T. Turki (2019) Tensor Decomposition-Based Unsupervised Feature Extraction Applied to Single-Cell Gene Expression Analysis. Frontiers in Genetics, 10(864): 10:3389/fgene.2019.00864
data(mouse_mid_brain)
data(mouse_mid_brain)
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(darr, ranks=NULL, max_iter=25, tol=1e-05) ## S4 method for signature 'DelayedArray' mpca(darr, ranks, max_iter, tol)
mpca(darr, ranks=NULL, max_iter=25, tol=1e-05) ## S4 method for signature 'DelayedArray' mpca(darr, ranks, max_iter, tol)
darr |
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 |
This function is an extension of the mpca
by DelayedArray.
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 darr
after compression
norm_percent
the percent of Frobenius norm explained by the approximation
fnorm_resid
the Frobenius norm of the error
fnorm(est-darr)
all_resids
vector containing the Frobenius norm of error for all the iterations
The length of ranks
must match darr@num_modes-1
.
H. Lu, K. Plataniotis, A. Venetsanopoulos, "Mpca: Multilinear principal component analysis of tensor objects". IEEE Trans. Neural networks, 2008.
library("DelayedRandomArray") darr <- RandomUnifArray(c(3,4,5)) mpca(darr, ranks=c(1,2))
library("DelayedRandomArray") darr <- RandomUnifArray(c(3,4,5)) mpca(darr, ranks=c(1,2))
Returns the outer product between two Tensors
outerProd(darr1, darr2) ## S4 method for signature 'DelayedArray,DelayedArray' outerProd(darr1, darr2)
outerProd(darr1, darr2) ## S4 method for signature 'DelayedArray,DelayedArray' outerProd(darr1, darr2)
darr1 |
first DelayedArray object |
darr2 |
second DelayedArray object |
NOTE: Sparse mode of outerProd is not available for now.
outer product between darr1
and darr2
library("DelayedRandomArray") darr1 <- RandomUnifArray(c(2,3)) darr2 <- RandomUnifArray(c(4,5)) outerProd(darr1, darr2)
library("DelayedRandomArray") darr1 <- RandomUnifArray(c(2,3)) darr2 <- RandomUnifArray(c(4,5)) outerProd(darr1, darr2)
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(darr, uranks=NULL, wranks=NULL, a=NULL, b=NULL) ## S4 method for signature 'DelayedArray' pvd(darr, uranks, wranks, a, b)
pvd(darr, uranks=NULL, wranks=NULL, a=NULL, b=NULL) ## S4 method for signature 'DelayedArray' pvd(darr, uranks, wranks, a, b)
darr |
3D DelayedArray (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 |
This function is an extension of the pvd
by DelayedArray.
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 darr
after compression
norm_percent
the percent of Frobenius norm explained by the approximation
fnorm_resid
the Frobenius norm of the error fnorm(est-darr)
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.
library("DelayedRandomArray") darr <- RandomUnifArray(c(3,4,5)) pvd(darr, uranks=rep(2,5), wranks=rep(3,5), a=2, b=3)
library("DelayedRandomArray") darr <- RandomUnifArray(c(3,4,5)) pvd(darr, uranks=rep(2,5), wranks=rep(3,5), a=2, b=3)
Returns the binded DelayedArray in row space.
rbind_list(L)
rbind_list(L)
L |
list of 2D DelayedArray |
This is a wrapper function to modebind_list
,
when the DelayedArrays are 2D.
2D DelayedArray object
The dimensions of row in each DelayedArray must match.
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) dlizt <- list( 'darr1' = RandomUnifArray(c(2,3)), 'darr2' = RandomUnifArray(c(2,3))) rbind_list(dlizt)
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) dlizt <- list( 'darr1' = RandomUnifArray(c(2,3)), 'darr2' = RandomUnifArray(c(2,3))) rbind_list(dlizt)
The inverse operation to rs_unfold
.
rs_fold(mat, m = NULL, modes = NULL) ## S4 method for signature 'DelayedArray' rs_fold(mat, m, modes)
rs_fold(mat, m = NULL, modes = NULL) ## S4 method for signature 'DelayedArray' rs_fold(mat, m, modes)
mat |
DelayedArray object (only 2D) |
m |
the mode corresponding to rs_unfold |
modes |
the original modes of the DelayedArray |
This function is an extension of the rs_fold
by DelayedArray.
This is a wrapper function to fold
.
DelayedArray (higher than 2D)
T. Kolda, B. Bader, "Tensor decomposition and applications". SIAM Applied Mathematics and Applications 2009.
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) matT2 <- rs_unfold(darr, m=2) identical( as.array(rs_fold(matT2, m=2, modes=c(2,3,4))), as.array(darr))
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) matT2 <- rs_unfold(darr, m=2) identical( as.array(rs_fold(matT2, m=2, modes=c(2,3,4))), as.array(darr))
Please see k_unfold
and unfold
.
rs_unfold(darr, m) ## S4 method for signature 'DelayedArray' rs_unfold(darr, m)
rs_unfold(darr, m) ## S4 method for signature 'DelayedArray' rs_unfold(darr, m)
darr |
DelayedArray object |
m |
mode to be unfolded on |
This function is an extension of the rs_unfold
by DelayedArray.
This is a wrapper function to unfold
.
See also rs_unfold(darr, m=NULL)
DelayedArray (2D)
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) matT2 <- rs_unfold(darr, m=2)
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) matT2 <- rs_unfold(darr, m=2)
Set whether the intermediate and output DelayedArray used in DelayedTensor is used as sparse tensor or not.
NOTE: Sparse mode is experimental! Whether it contributes to higher speed and lower memory is quite dependent on the sparsity of the DelayedArray, and the current implementation does not recognize the block size, which may cause Out-of-Memory errors.
setSparse(as.sparse=FALSE)
setSparse(as.sparse=FALSE)
as.sparse |
TRUE or FALSE (Default: FALSE) |
Nothing
setSparse(TRUE) setSparse(FALSE)
setSparse(TRUE) setSparse(FALSE)
Set the verbose message to monitor the block-processing procedure.
setVerbose(as.verbose=FALSE)
setVerbose(as.verbose=FALSE)
as.verbose |
TRUE or FALSE (Default: FALSE) |
Nothing
setVerbose(TRUE) setVerbose(FALSE)
setVerbose(TRUE) setVerbose(FALSE)
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(darr, list_mat, ms=NULL)
ttl(darr, list_mat, ms=NULL)
darr |
DelayedArray object with K modes |
list_mat |
a list of 2D DelayedArray objects |
ms |
a vector of modes to contract on
(order should match the order of |
This function is an extension of the ttl
by DelayedArray.
This is a wrapper function to unfold
.
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 darr
with
three matrices mat1
, mat2
, mat3
on modes 1, 2, and 3.
We could do ttm(ttm(ttm(darr,mat1,1),mat2,2),3)
,
or we could do ttl(darr,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).
DelayedArray object with K modes (Tensor)
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.
library("DelayedRandomArray") darr <- RandomUnifArray(c(3,4,5)) dlizt <- list( 'darr1' = RandomUnifArray(c(10,3)), 'darr2' = RandomUnifArray(c(10,4))) ttl(darr, dlizt, ms=c(1,2))
library("DelayedRandomArray") darr <- RandomUnifArray(c(3,4,5)) dlizt <- list( 'darr1' = RandomUnifArray(c(10,3)), 'darr2' = RandomUnifArray(c(10,4))) ttl(darr, dlizt, ms=c(1,2))
Contracted (m-Mode) product between a DelayedArray (Tensor) of arbitrary number of modes and a matrix. The result is folded back into Tensor.
ttm(darr, mat, m = NULL) ## S4 method for signature 'DelayedArray,DelayedArray' ttm(darr, mat, m)
ttm(darr, mat, m = NULL) ## S4 method for signature 'DelayedArray,DelayedArray' ttm(darr, mat, m)
darr |
DelayedArray object |
mat |
input 2D DelayedArray with same number columns as
the |
m |
the mode to contract on |
This function is an extension of the ttm
by DelayedArray.
By definition, rs_unfold(ttm(darr, mat), m) = mat%*%rs_unfold(darr, m)
,
so the number of columns in mat
must match the
m
th mode of darr
. For the math on the m-Mode Product,
see Kolda and Bader (2009).
a DelayedArray object with K modes
The m
th mode of darr
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.
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) mat <- RandomUnifArray(c(10,4)) ttm(darr, mat, m=3)
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) mat <- RandomUnifArray(c(10,4)) ttm(darr, 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 = darr@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(darr, ranks=NULL, max_iter=25, tol=1e-05) ## S4 method for signature 'DelayedArray' tucker(darr, ranks, max_iter, tol)
tucker(darr, ranks=NULL, max_iter=25, tol=1e-05) ## S4 method for signature 'DelayedArray' tucker(darr, ranks, max_iter, tol)
darr |
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 |
This function is an extension of the tucker
by DelayedArray.
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 darr
after compression
norm_percent
the percent of Frobenius norm explained by the approximation
fnorm_resid
the Frobenius norm of the error
fnorm(est-darr)
all_resids
vector containing the Frobenius norm of error for all the iterations
The length of ranks
must match darr@num_modes
.
T. Kolda, B. Bader, "Tensor decomposition and applications". SIAM Applied Mathematics and Applications 2009.
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) tucker(darr, ranks=c(1,2,3))
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) tucker(darr, ranks=c(1,2,3))
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 dim(darr)
.
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(darr, row_idx, col_idx) ## S4 method for signature 'DelayedArray' unfold(darr, row_idx, col_idx)
unfold(darr, row_idx, col_idx) ## S4 method for signature 'DelayedArray' unfold(darr, row_idx, col_idx)
darr |
DelayedArray object |
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 |
This function is an extension of the unfold
by DelayedArray.
For Row Space Unfolding or m-mode Unfolding,
see rs_unfold
.
For Column Space Unfolding or matvec,
see cs_unfold
.
vec
returns the vectorization of the tensor.
2D DelayedArray 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
, matvec
,
rs_unfold
, cs_unfold
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) unfold(darr, row_idx=2, col_idx=c(3,1))
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) unfold(darr, 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) ## S4 method for signature 'DelayedArray' unmatvec(mat, modes)
unmatvec(mat, modes = NULL) ## S4 method for signature 'DelayedArray' unmatvec(mat, modes)
mat |
DelayedArray object (only 2D) |
modes |
the modes of the output DelayedArray |
This function is an extension of the unmatvec
by DelayedArray.
This is a wrapper function to fold
.
DelayedArray object with modes given by modes
T. Kolda, B. Bader, "Tensor decomposition and applications". SIAM Applied Mathematics and Applications 2009.
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) matT1 <- matvec(darr) identical( as.array(unmatvec(matT1, modes=c(2,3,4))), as.array(darr))
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) matT1 <- matvec(darr) identical( as.array(unmatvec(matT1, modes=c(2,3,4))), as.array(darr))
Change the dimension of DelayedArray from multi-dimension (e.g. array) to single-dimension (e.g. vector).
vec(darr) ## S4 method for signature 'DelayedArray' vec(darr)
vec(darr) ## S4 method for signature 'DelayedArray' vec(darr)
darr |
DelayedArray object |
This function is an extension of the vec
by DelayedArray.
1D DelayedArray (vector) with length prod(dim(darr))
T. Kolda, B. Bader, "Tensor decomposition and applications". SIAM Applied Mathematics and Applications 2009.
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) vec(darr)
library("DelayedRandomArray") darr <- RandomUnifArray(c(2,3,4)) vec(darr)