DelayedTensor 1.4.0
Authors: Koki Tsuyuzaki [aut, cre]
Last modified: 2022-11-01 15:19:15
Compiled: Tue Nov 1 17:20:04 2022
suppressPackageStartupMessages(library("DelayedTensor"))
suppressPackageStartupMessages(library("DelayedArray"))
suppressPackageStartupMessages(library("HDF5Array"))
suppressPackageStartupMessages(library("DelayedRandomArray"))
darr1 <- RandomUnifArray(c(2,3,4))
darr2 <- RandomUnifArray(c(2,3,4))
There are several settings in DelayedTensor.
First, the sparsity of the intermediate DelayedArray objects
calculated inside DelayedTensor is set by setSparse
.
Note that the 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, when the data is extremely huge.
Here, we specify as.sparse
as FALSE
(this is also the default value for now).
DelayedTensor::setSparse(as.sparse=FALSE)
Next, the verbose message is suppressed by setVerbose
.
This is useful when we want to monitor the calculation process.
Here we specify as.verbose
as FALSE
(this is also the default value for now).
DelayedTensor::setVerbose(as.verbose=FALSE)
The block size of block processing is specified by setAutoBlockSize
.
When the sparse mode is off, all the functions of DelayedTensor
are performed as block processing,
in which each block vector/matrix/tensor is expanded to memory space
from on-disk file incrementally so as not to exceed the specified size.
Here, we specify the block size as 1E+8
.
setAutoBlockSize(size=1E+8)
## automatic block size set to 1e+08 bytes (was 1e+08)
Finally, the temporal directory to store the intermediate HDF5 files during running DelayedTensor is specified by setHDF5DumpDir
.
Note that in many systems the /var
directory has the storage limitation, so if there is no enough space, user should specify the other directory.
# tmpdir <- paste(sample(c(letters,1:9), 10), collapse="")
# dir.create(tmpdir, recursive=TRUE))
tmpdir <- tempdir()
setHDF5DumpDir(tmpdir)
These specified values are also extracted by each getter function.
DelayedTensor::getSparse()
## $delayedtensor.sparse
## [1] FALSE
DelayedTensor::getVerbose()
## $delayedtensor.verbose
## [1] FALSE
getAutoBlockSize()
## [1] 1e+08
getHDF5DumpDir()
## [1] "/tmp/RtmpGvxYlP"
Unfold (a.k.a. matricizing) operations are used to reshape a tensor into a matrix.
In unfold
, row_idx
and col_idx
are specified to set which modes are used
as the row/column.
dmat1 <- DelayedTensor::unfold(darr1, row_idx=c(1,2), col_idx=3)
dmat1
## <6 x 4> matrix of class HDF5Matrix and type "double":
## [,1] [,2] [,3] [,4]
## [1,] 0.06169532 0.87277098 0.57033972 0.42749688
## [2,] 0.84892416 0.13521219 0.35014967 0.72173113
## [3,] 0.78152225 0.07971139 0.39378891 0.63189357
## [4,] 0.79607687 0.99588605 0.69279147 0.16429066
## [5,] 0.46932250 0.35048294 0.74312051 0.68148783
## [6,] 0.72767438 0.28030408 0.60529128 0.03640777
fold
is the inverse operation of unfold
, which is used to reshape
a matrix into a tensor.
In fold
, row_idx
/col_idx
are specified to set which modes correspond
the row/column of the output tensor and modes
is specified to set the mode of the output tensor.
dmat1_to_darr1 <- DelayedTensor::fold(dmat1,
row_idx=c(1,2), col_idx=3, modes=dim(darr1))
dmat1_to_darr1
## <2 x 3 x 4> array of class DelayedArray and type "double":
## ,,1
## [,1] [,2] [,3]
## [1,] 0.06169532 0.78152225 0.46932250
## [2,] 0.84892416 0.79607687 0.72767438
##
## ,,2
## [,1] [,2] [,3]
## [1,] 0.87277098 0.07971139 0.35048294
## [2,] 0.13521219 0.99588605 0.28030408
##
## ,,3
## [,1] [,2] [,3]
## [1,] 0.5703397 0.3937889 0.7431205
## [2,] 0.3501497 0.6927915 0.6052913
##
## ,,4
## [,1] [,2] [,3]
## [1,] 0.42749688 0.63189357 0.68148783
## [2,] 0.72173113 0.16429066 0.03640777
identical(as.array(darr1), as.array(dmat1_to_darr1))
## [1] TRUE
There are some wrapper functions of unfold
and fold
.
For example, in k_unfold
, mode m
is used as the row, and the other modes
are is used as the column.
k_fold
is the inverse operation of k_unfold
.
dmat2 <- DelayedTensor::k_unfold(darr1, m=1)
dmat2_to_darr1 <- k_fold(dmat2, m=1, modes=dim(darr1))
identical(as.array(darr1), as.array(dmat2_to_darr1))
## [1] TRUE
dmat3 <- DelayedTensor::k_unfold(darr1, m=2)
dmat3_to_darr1 <- k_fold(dmat3, m=2, modes=dim(darr1))
identical(as.array(darr1), as.array(dmat3_to_darr1))
## [1] TRUE
dmat4 <- DelayedTensor::k_unfold(darr1, m=3)
dmat4_to_darr1 <- k_fold(dmat4, m=3, modes=dim(darr1))
identical(as.array(darr1), as.array(dmat4_to_darr1))
## [1] TRUE
In rs_unfold
, mode m
is used as the row, and the other modes
are is used as the column.
rs_fold
and rs_unfold
also perform the same operations.
On the other hand, cs_unfold
specifies the mode m
as the column
and the other modes are specified as the column.
cs_fold
is the inverse operation of cs_unfold
.
dmat8 <- DelayedTensor::cs_unfold(darr1, m=1)
dmat8_to_darr1 <- DelayedTensor::cs_fold(dmat8, m=1, modes=dim(darr1))
identical(as.array(darr1), as.array(dmat8_to_darr1))
## [1] TRUE
dmat9 <- DelayedTensor::cs_unfold(darr1, m=2)
dmat9_to_darr1 <- DelayedTensor::cs_fold(dmat9, m=2, modes=dim(darr1))
identical(as.array(darr1), as.array(dmat9_to_darr1))
## [1] TRUE
dmat10 <- DelayedTensor::cs_unfold(darr1, m=3)
dmat10_to_darr1 <- DelayedTensor::cs_fold(dmat10, m=3, modes=dim(darr1))
identical(as.array(darr1), as.array(dmat10_to_darr1))
## [1] TRUE
In matvec
, m=2 is specified as unfold.
unmatvec
is the inverse operation of matvec
.
dmat11 <- DelayedTensor::matvec(darr1)
dmat11_darr1 <- DelayedTensor::unmatvec(dmat11, modes=dim(darr1))
identical(as.array(darr1), as.array(dmat11_darr1))
## [1] TRUE
ttm
multiplies a tensor by a matrix.
m
specifies in which mode the matrix will be multiplied.
dmatZ <- RandomUnifArray(c(10,4))
DelayedTensor::ttm(darr1, dmatZ, m=3)
## <2 x 3 x 10> array of class DelayedArray and type "double":
## ,,1
## [,1] [,2] [,3]
## [1,] 1.523594 1.243578 1.626893
## [2,] 1.360216 1.839284 1.054345
##
## ,,2
## [,1] [,2] [,3]
## [1,] 1.2525427 1.0133236 1.1667560
## [2,] 1.1707036 1.5856715 0.7412799
##
## ,,3
## [,1] [,2] [,3]
## [1,] 0.8364524 0.7030880 0.8322550
## [2,] 0.7763512 1.2651599 0.7384567
##
## ...
##
## ,,8
## [,1] [,2] [,3]
## [1,] 0.7397048 0.9649453 1.0112964
## [2,] 1.0786792 0.8882292 0.5295652
##
## ,,9
## [,1] [,2] [,3]
## [1,] 0.7592775 0.4278659 0.7142002
## [2,] 0.4552795 0.8320064 0.4573291
##
## ,,10
## [,1] [,2] [,3]
## [1,] 0.3000667 0.3805549 0.4665307
## [2,] 0.3814704 0.5307788 0.4524404
ttl
multiplies a tensor by multiple matrices.
ms
specifies in which mode these matrices will be multiplied.
dmatX <- RandomUnifArray(c(10,2))
dmatY <- RandomUnifArray(c(10,3))
dlizt <- list(dmatX = dmatX, dmatY = dmatY)
DelayedTensor::ttl(darr1, dlizt, ms=c(1,2))
## <10 x 10 x 4> array of class DelayedArray and type "double":
## ,,1
## [,1] [,2] [,3] ... [,9] [,10]
## [1,] 0.8138620 1.1081015 0.9805884 . 1.147542 1.678992
## [2,] 0.8365134 1.1211665 0.9956324 . 1.150301 1.683906
## ... . . . . . .
## [9,] 0.6433780 0.8904178 0.7851261 . 0.9308568 1.3612449
## [10,] 0.9827818 1.3716736 1.2072513 . 1.4408419 2.1064715
##
## ...
##
## ,,4
## [,1] [,2] [,3] ... [,9] [,10]
## [1,] 0.7837297 0.9383808 0.8661579 . 0.8535495 1.2700621
## [2,] 0.7560630 0.9021131 0.8393315 . 0.7956859 1.1931488
## ... . . . . . .
## [9,] 0.6597419 0.7924783 0.7260841 . 0.7410385 1.0951877
## [10,] 1.0398730 1.2510032 1.1421559 . 1.1849038 1.7457546
vec
collapses a DelayedArray into
a 1D DelayedArray (vector).
DelayedTensor::vec(darr1)
## <24> array of class HDF5Array and type "double":
## [1] [2] [3] . [23] [24]
## 0.06169532 0.84892416 0.78152225 . 0.68148783 0.03640777
fnorm
calculates the Frobenius norm of a DelayedArray.
DelayedTensor::fnorm(darr1)
## [1] 2.880198
innerProd
calculates the inner product value of two
DelayedArray.
DelayedTensor::innerProd(darr1, darr2)
## [1] 5.428505
Inner product multiplies two tensors and collapses to 0D tensor (norm). On the other hand, the outer product is an operation that leaves all subscripts intact.
DelayedTensor::outerProd(darr1[,,1], darr2[,,1])
## <2 x 3 x 2 x 3> array of class HDF5Array and type "double":
## ,,1,1
## [,1] [,2] [,3]
## [1,] 0.01594678 0.20200504 0.12130878
## [2,] 0.21942685 0.20576707 0.18808664
##
## ,,2,1
## [,1] [,2] [,3]
## [1,] 0.01410503 0.17867476 0.10729840
## [2,] 0.19408446 0.18200230 0.16636385
##
## ,,1,2
## [,1] [,2] [,3]
## [1,] 0.01051657 0.13321806 0.08000058
## [2,] 0.14470737 0.13569903 0.12403917
##
## ,,2,2
## [,1] [,2] [,3]
## [1,] 0.03428019 0.43424250 0.26077284
## [2,] 0.47169348 0.44232958 0.40432265
##
## ,,1,3
## [,1] [,2] [,3]
## [1,] 0.05819237 0.73714874 0.44267517
## [2,] 0.80072370 0.75087698 0.68635827
##
## ,,2,3
## [,1] [,2] [,3]
## [1,] 0.03862502 0.48928042 0.29382441
## [2,] 0.53147812 0.49839250 0.45556839
Using DelayedDiagonalArray
, we can originally create a diagonal
DelayedArray by specifying the dimensions (modes) and the values.
dgdarr <- DelayedTensor::DelayedDiagonalArray(c(5,6,7), 1:5)
dgdarr
## <5 x 6 x 7> sparse array of class DelayedArray and type "integer":
## ,,1
## [,1] [,2] [,3] [,4] [,5] [,6]
## [1,] 1 0 0 0 0 0
## [2,] 0 0 0 0 0 0
## [3,] 0 0 0 0 0 0
## [4,] 0 0 0 0 0 0
## [5,] 0 0 0 0 0 0
##
## ...
##
## ,,7
## [,1] [,2] [,3] [,4] [,5] [,6]
## [1,] 0 0 0 0 0 0
## [2,] 0 0 0 0 0 0
## [3,] 0 0 0 0 0 0
## [4,] 0 0 0 0 0 0
## [5,] 0 0 0 0 0 0
Similar to the diag
of the base package,
the diag
of DelayedTensor is used to extract
and assign values to DelayedArray.
DelayedTensor::diag(dgdarr)
## <5> array of class DelayedArray and type "integer":
## [1] [2] [3] [4] [5]
## 1 2 3 4 5
DelayedTensor::diag(dgdarr) <- c(1111, 2222, 3333, 4444, 5555)
DelayedTensor::diag(dgdarr)
## <5> array of class DelayedArray and type "double":
## [1] [2] [3] [4] [5]
## 1111 2222 3333 4444 5555
modeSum
calculates the summation for a given mode m
of
a DelayedArray.
The mode specified as m
is collapsed into 1D as follows.
DelayedTensor::modeSum(darr1, m=1)
## <1 x 3 x 4> array of class DelayedArray and type "double":
## ,,1
## [,1] [,2] [,3]
## [1,] 0.9106195 1.5775991 1.1969969
##
## ,,2
## [,1] [,2] [,3]
## [1,] 1.007983 1.075597 0.630787
##
## ,,3
## [,1] [,2] [,3]
## [1,] 0.9204894 1.0865804 1.3484118
##
## ,,4
## [,1] [,2] [,3]
## [1,] 1.1492280 0.7961842 0.7178956
DelayedTensor::modeSum(darr1, m=2)
## <2 x 1 x 4> array of class DelayedArray and type "double":
## ,,1
## [,1]
## [1,] 1.312540
## [2,] 2.372675
##
## ,,2
## [,1]
## [1,] 1.302965
## [2,] 1.411402
##
## ,,3
## [,1]
## [1,] 1.707249
## [2,] 1.648232
##
## ,,4
## [,1]
## [1,] 1.7408783
## [2,] 0.9224296
DelayedTensor::modeSum(darr1, m=3)
## <2 x 3 x 1> array of class DelayedArray and type "double":
## ,,1
## [,1] [,2] [,3]
## [1,] 1.932303 1.886916 2.244414
## [2,] 2.056017 2.649045 1.649678
Similar to modeSum
, modeMean
calculates the average value
for a given mode m
of a DelayedArray.
DelayedTensor::modeMean(darr1, m=1)
## <1 x 3 x 4> array of class DelayedArray and type "double":
## ,,1
## [,1] [,2] [,3]
## [1,] 0.4553097 0.7887996 0.5984984
##
## ,,2
## [,1] [,2] [,3]
## [1,] 0.5039916 0.5377987 0.3153935
##
## ,,3
## [,1] [,2] [,3]
## [1,] 0.4602447 0.5432902 0.6742059
##
## ,,4
## [,1] [,2] [,3]
## [1,] 0.5746140 0.3980921 0.3589478
DelayedTensor::modeMean(darr1, m=2)
## <2 x 1 x 4> array of class DelayedArray and type "double":
## ,,1
## [,1]
## [1,] 0.4375134
## [2,] 0.7908918
##
## ,,2
## [,1]
## [1,] 0.4343218
## [2,] 0.4704674
##
## ,,3
## [,1]
## [1,] 0.5690830
## [2,] 0.5494108
##
## ,,4
## [,1]
## [1,] 0.5802928
## [2,] 0.3074765
DelayedTensor::modeMean(darr1, m=3)
## <2 x 3 x 1> array of class DelayedArray and type "double":
## ,,1
## [,1] [,2] [,3]
## [1,] 0.4830757 0.4717290 0.5611034
## [2,] 0.5140043 0.6622613 0.4124194
There are some tensor specific product such as Hadamard product, Kronecker product, and Khatri-Rao product.
Suppose a tensor \(A \in \Re ^{I \times J}\) and a tensor \(B \in \Re ^{I \times J}\).
Hadamard product is defined as the element-wise product of \(A\) and \(B\).
Hadamard product can be extended to higher-order tensors.
\[ A \circ B = \begin{bmatrix} a_{11}b_{11} & a_{12}b_{12} & \cdots & a_{1J}b_{1J} \\ a_{21}b_{21} & a_{22}b_{22} & \cdots & a_{2J}b_{2J} \\ \vdots & \vdots & \ddots & \vdots \\ a_{I1}b_{I1} & a_{I2}b_{I2} & \cdots & a_{IJ}b_{IJ} \\ \end{bmatrix} \]
hadamard
calculates Hadamard product of two DelayedArray
objects.
prod_h <- DelayedTensor::hadamard(darr1, darr2)
dim(prod_h)
## [1] 2 3 4
hadamard_list
calculates Hadamard product of multiple
DelayedArray objects.
prod_hl <- DelayedTensor::hadamard_list(list(darr1, darr2))
dim(prod_hl)
## [1] 2 3 4
Suppose a tensor \(A \in \Re ^{I \times J}\) and a tensor \(B \in \Re ^{K \times L}\).
Kronecker product is defined as all the possible combination of element-wise product and the dimensions of output tensor are \({IK \times JL}\).
Kronecker product can be extended to higher-order tensors.
\[ A \otimes B = \begin{bmatrix} a_{11}B & a_{12}B & \cdots & a_{1J}B \\ a_{21}B & a_{22}B & \cdots & a_{2J}B \\ \vdots & \vdots & \ddots & \vdots \\ a_{I1}B & a_{I2}B & \cdots & a_{IJ}B \\ \end{bmatrix} \]
kronecker
calculates Kronecker product of two DelayedArray
objects.
prod_kron <- DelayedTensor::kronecker(darr1, darr2)
dim(prod_kron)
## [1] 4 9 16
kronecker_list
calculates Kronecker product of multiple
DelayedArray objects.
prod_kronl <- DelayedTensor::kronecker_list(list(darr1, darr2))
dim(prod_kronl)
## [1] 4 9 16
Suppose a tensor \(A \in \Re ^{I \times J}\) and a tensor \(B \in \Re ^{K \times J}\).
Khatri-Rao product is defined as the column-wise Kronecker product and the dimensions of output tensor is \({IK \times J}\).
\[ A \odot B = \begin{bmatrix} a_{1} \otimes a_{1} & a_{2} \otimes a_{2} & \cdots & a_{J} \otimes a_{J} \\ \end{bmatrix} \]
Khatri-Rao product can only be used for 2D tensors (matrices).
khatri_rao
calculates Khatri-Rao product of two DelayedArray
objects.
prod_kr <- DelayedTensor::khatri_rao(darr1[,,1], darr2[,,1])
dim(prod_kr)
## [1] 4 3
khatri_rao_list
calculates Khatri-Rao product of multiple
DelayedArray objects.
prod_krl <- DelayedTensor::khatri_rao_list(list(darr1[,,1], darr2[,,1]))
dim(prod_krl)
## [1] 4 3
list_rep
replicates an arbitrary number of any R object.
str(DelayedTensor::list_rep(darr1, 3))
## List of 3
## $ :Formal class 'RandomUnifArray' [package "DelayedRandomArray"] with 1 slot
## .. ..@ seed:Formal class 'RandomUnifArraySeed' [package "DelayedRandomArray"] with 6 slots
## .. .. .. ..@ min : num 0
## .. .. .. ..@ max : num 1
## .. .. .. ..@ dim : int [1:3] 2 3 4
## .. .. .. ..@ chunkdim: int [1:3] 2 3 4
## .. .. .. ..@ seeds :List of 1
## .. .. .. .. ..$ : int [1:2] 846547912 495635619
## .. .. .. ..@ sparse : logi FALSE
## $ :Formal class 'RandomUnifArray' [package "DelayedRandomArray"] with 1 slot
## .. ..@ seed:Formal class 'RandomUnifArraySeed' [package "DelayedRandomArray"] with 6 slots
## .. .. .. ..@ min : num 0
## .. .. .. ..@ max : num 1
## .. .. .. ..@ dim : int [1:3] 2 3 4
## .. .. .. ..@ chunkdim: int [1:3] 2 3 4
## .. .. .. ..@ seeds :List of 1
## .. .. .. .. ..$ : int [1:2] 846547912 495635619
## .. .. .. ..@ sparse : logi FALSE
## $ :Formal class 'RandomUnifArray' [package "DelayedRandomArray"] with 1 slot
## .. ..@ seed:Formal class 'RandomUnifArraySeed' [package "DelayedRandomArray"] with 6 slots
## .. .. .. ..@ min : num 0
## .. .. .. ..@ max : num 1
## .. .. .. ..@ dim : int [1:3] 2 3 4
## .. .. .. ..@ chunkdim: int [1:3] 2 3 4
## .. .. .. ..@ seeds :List of 1
## .. .. .. .. ..$ : int [1:2] 846547912 495635619
## .. .. .. ..@ sparse : logi FALSE
modebind_list
collapses multiple DelayedArray objects
into single DelayedArray object.
m
specifies the collapsed dimension.
dim(DelayedTensor::modebind_list(list(darr1, darr2), m=1))
## [1] 4 3 4
dim(DelayedTensor::modebind_list(list(darr1, darr2), m=2))
## [1] 2 6 4
dim(DelayedTensor::modebind_list(list(darr1, darr2), m=3))
## [1] 2 3 8
rbind_list
is the row-wise modebind_list
and
collapses multiple 2D DelayedArray objects
into single DelayedArray object.
dim(DelayedTensor::rbind_list(list(darr1[,,1], darr2[,,1])))
## [1] 4 3
cbind_list
is the column-wise modebind_list
and
collapses multiple 2D DelayedArray objects
into single DelayedArray object.
dim(DelayedTensor::cbind_list(list(darr1[,,1], darr2[,,1])))
## [1] 2 6
## R version 4.2.1 (2022-06-23)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.5 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.16-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.16-bioc/R/lib/libRlapack.so
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_GB LC_COLLATE=C
## [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
## [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
## [9] LC_ADDRESS=C LC_TELEPHONE=C
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## attached base packages:
## [1] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] DelayedRandomArray_1.6.0 HDF5Array_1.26.0 rhdf5_2.42.0
## [4] DelayedArray_0.24.0 IRanges_2.32.0 S4Vectors_0.36.0
## [7] MatrixGenerics_1.10.0 matrixStats_0.62.0 BiocGenerics_0.44.0
## [10] Matrix_1.5-1 DelayedTensor_1.4.0 BiocStyle_2.26.0
##
## loaded via a namespace (and not attached):
## [1] Rcpp_1.0.9 rTensor_1.4.8 bslib_0.4.0
## [4] compiler_4.2.1 BiocManager_1.30.19 jquerylib_0.1.4
## [7] rhdf5filters_1.10.0 tools_4.2.1 digest_0.6.30
## [10] jsonlite_1.8.3 evaluate_0.17 lattice_0.20-45
## [13] rlang_1.0.6 cli_3.4.1 parallel_4.2.1
## [16] yaml_2.3.6 xfun_0.34 fastmap_1.1.0
## [19] stringr_1.4.1 knitr_1.40 sass_0.4.2
## [22] grid_4.2.1 R6_2.5.1 BiocParallel_1.32.0
## [25] rmarkdown_2.17 bookdown_0.29 irlba_2.3.5.1
## [28] BiocSingular_1.14.0 Rhdf5lib_1.20.0 magrittr_2.0.3
## [31] codetools_0.2-18 htmltools_0.5.3 rsvd_1.0.5
## [34] beachmat_2.14.0 dqrng_0.3.0 ScaledMatrix_1.6.0
## [37] stringi_1.7.8 einsum_0.1.0 cachem_1.0.6