Authors: Koki Tsuyuzaki [aut, cre]
Last modified: 2024-04-30 20:49:30.349857
Compiled: Tue Apr 30 23:25:38 2024

1 Setting

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/RtmpdkJd2z"

2 Tensor Arithmetic Operations

2.1 Unfold/Fold Operations

Unfold (a.k.a. matricizing) operations are used to reshape a tensor into a matrix.

Figure 1: Unfold/Fold Operasions

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> HDF5Matrix object of type "double":
##            [,1]       [,2]       [,3]       [,4]
## [1,] 0.43021098 0.06662193 0.07909128 0.32360773
## [2,] 0.17419051 0.05271241 0.69457602 0.23985783
## [3,] 0.53236558 0.23141741 0.87853694 0.02254453
## [4,] 0.47042131 0.43398851 0.53227820 0.87985983
## [5,] 0.70151454 0.33933496 0.51644139 0.62438232
## [6,] 0.10665020 0.31976206 0.19971189 0.24892296

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> DelayedArray object of type "double":
## ,,1
##           [,1]      [,2]      [,3]
## [1,] 0.4302110 0.5323656 0.7015145
## [2,] 0.1741905 0.4704213 0.1066502
## 
## ,,2
##            [,1]       [,2]       [,3]
## [1,] 0.06662193 0.23141741 0.33933496
## [2,] 0.05271241 0.43398851 0.31976206
## 
## ,,3
##            [,1]       [,2]       [,3]
## [1,] 0.07909128 0.87853694 0.51644139
## [2,] 0.69457602 0.53227820 0.19971189
## 
## ,,4
##            [,1]       [,2]       [,3]
## [1,] 0.32360773 0.02254453 0.62438232
## [2,] 0.23985783 0.87985983 0.24892296
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> DelayedArray object of type "double":
## ,,1
##           [,1]      [,2]      [,3]
## [1,] 0.6331984 0.7744297 1.2612956
## [2,] 0.5232865 1.2653152 0.3943736
## 
## ,,2
##           [,1]      [,2]      [,3]
## [1,] 0.3326510 0.3993356 0.5634862
## [2,] 0.1473689 0.4205377 0.1108611
## 
## ,,3
##           [,1]      [,2]      [,3]
## [1,] 0.4675453 0.6149040 1.0940616
## [2,] 0.6293886 1.3249467 0.4310947
## 
## ...
## 
## ,,8
##           [,1]      [,2]      [,3]
## [1,] 0.5869515 1.1764349 1.4906429
## [2,] 0.8019484 1.5953386 0.6369104
## 
## ,,9
##           [,1]      [,2]      [,3]
## [1,] 0.4887139 0.8737401 1.2810867
## [2,] 0.7723742 1.5307141 0.5683950
## 
## ,,10
##           [,1]      [,2]      [,3]
## [1,] 0.2759367 0.3207960 0.5948424
## [2,] 0.1645167 0.6301346 0.2654136

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> DelayedArray object of type "double":
## ,,1
##            [,1]      [,2]      [,3] ...      [,9]     [,10]
##  [1,] 0.4100945 0.8605137 0.8395522   . 0.8600078 0.3930065
##  [2,] 0.2743268 0.7301880 0.6255925   . 0.6913216 0.3362309
##   ...         .         .         .   .         .         .
##  [9,] 0.3548820 0.8517110 0.7708384   . 0.8245880 0.3908879
## [10,] 0.3093642 0.6831756 0.6474224   . 0.6743118 0.3126186
## 
## ...
## 
## ,,4
##            [,1]      [,2]      [,3] ...      [,9]     [,10]
##  [1,] 0.4020881 0.7848356 0.8053501   . 0.8101782 0.3589531
##  [2,] 0.3607688 1.0313360 0.8790178   . 1.0082173 0.4822665
##   ...         .         .         .   .         .         .
##  [9,] 0.4115346 1.0298660 0.9326151   . 1.0240424 0.4783440
## [10,] 0.3235343 0.7035310 0.6824515   . 0.7138268 0.3240954

2.2 Vectorization

vec collapses a DelayedArray into a 1D DelayedArray (vector).

Figure 2: Vectorization

DelayedTensor::vec(darr1)
## <24> HDF5Array object of type "double":
##       [1]       [2]       [3]         .      [23]      [24] 
## 0.4302110 0.1741905 0.5323656         . 0.6243823 0.2489230

2.3 Norm Operations

fnorm calculates the Frobenius norm of a DelayedArray.

Figure 3: Norm Operations

DelayedTensor::fnorm(darr1)
## [1] 2.220205

innerProd calculates the inner product value of two DelayedArray.

DelayedTensor::innerProd(darr1, darr2)
## [1] 5.327103

2.4 Outer Product

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.

Figure 4: Outer Product

DelayedTensor::outerProd(darr1[,,1], darr2[,,1])
## <2 x 3 x 2 x 3> HDF5Array object of type "double":
## ,,1,1
##           [,1]      [,2]      [,3]
## [1,] 0.3479442 0.4305644 0.5673680
## [2,] 0.1408811 0.3804654 0.0862561
## 
## ,,2,1
##            [,1]       [,2]       [,3]
## [1,] 0.27109551 0.33546777 0.44205622
## [2,] 0.10976536 0.29643386 0.06720514
## 
## ,,1,2
##            [,1]       [,2]       [,3]
## [1,] 0.09213031 0.11400686 0.15023036
## [2,] 0.03730315 0.10074141 0.02283930
## 
## ,,2,2
##            [,1]       [,2]       [,3]
## [1,] 0.22752650 0.28155319 0.37101132
## [2,] 0.09212446 0.24879261 0.05640429
## 
## ,,1,3
##            [,1]       [,2]       [,3]
## [1,] 0.19950535 0.24687836 0.32531922
## [2,] 0.08077883 0.21815242 0.04945779
## 
## ,,2,3
##            [,1]       [,2]       [,3]
## [1,] 0.05713668 0.07070391 0.09316873
## [2,] 0.02313439 0.06247704 0.01416430

2.5 Diagonal Operations

Using DelayedDiagonalArray, we can originally create a diagonal DelayedArray by specifying the dimensions (modes) and the values.

Figure 5: Diagonal Operations

dgdarr <- DelayedTensor::DelayedDiagonalArray(c(5,6,7), 1:5)
dgdarr
## <5 x 6 x 7> sparse DelayedArray object of 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> DelayedArray object of 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> DelayedArray object of type "double":
##  [1]  [2]  [3]  [4]  [5] 
## 1111 2222 3333 4444 5555

2.6 Mode-wise Operations

modeSum calculates the summation for a given mode m of a DelayedArray. The mode specified as m is collapsed into 1D as follows.

Figure 6: Mode-wise Operations

DelayedTensor::modeSum(darr1, m=1)
## <1 x 3 x 4> DelayedArray object of type "double":
## ,,1
##           [,1]      [,2]      [,3]
## [1,] 0.6044015 1.0027869 0.8081647
## 
## ,,2
##           [,1]      [,2]      [,3]
## [1,] 0.1193343 0.6654059 0.6590970
## 
## ,,3
##           [,1]      [,2]      [,3]
## [1,] 0.7736673 1.4108151 0.7161533
## 
## ,,4
##           [,1]      [,2]      [,3]
## [1,] 0.5634656 0.9024044 0.8733053
DelayedTensor::modeSum(darr1, m=2)
## <2 x 1 x 4> DelayedArray object of type "double":
## ,,1
##          [,1]
## [1,] 1.664091
## [2,] 0.751262
## 
## ,,2
##           [,1]
## [1,] 0.6373743
## [2,] 0.8064630
## 
## ,,3
##          [,1]
## [1,] 1.474070
## [2,] 1.426566
## 
## ,,4
##           [,1]
## [1,] 0.9705346
## [2,] 1.3686406
DelayedTensor::modeSum(darr1, m=3)
## <2 x 3 x 1> DelayedArray object of type "double":
## ,,1
##           [,1]      [,2]      [,3]
## [1,] 0.8995319 1.6648645 2.1816732
## [2,] 1.1613368 2.3165478 0.8750471

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> DelayedArray object of type "double":
## ,,1
##           [,1]      [,2]      [,3]
## [1,] 0.3022007 0.5013934 0.4040824
## 
## ,,2
##            [,1]       [,2]       [,3]
## [1,] 0.05966717 0.33270296 0.32954851
## 
## ,,3
##           [,1]      [,2]      [,3]
## [1,] 0.3868336 0.7054076 0.3580766
## 
## ,,4
##           [,1]      [,2]      [,3]
## [1,] 0.2817328 0.4512022 0.4366526
DelayedTensor::modeMean(darr1, m=2)
## <2 x 1 x 4> DelayedArray object of type "double":
## ,,1
##           [,1]
## [1,] 0.5546970
## [2,] 0.2504207
## 
## ,,2
##           [,1]
## [1,] 0.2124581
## [2,] 0.2688210
## 
## ,,3
##           [,1]
## [1,] 0.4913565
## [2,] 0.4755220
## 
## ,,4
##           [,1]
## [1,] 0.3235115
## [2,] 0.4562135
DelayedTensor::modeMean(darr1, m=3)
## <2 x 3 x 1> DelayedArray object of type "double":
## ,,1
##           [,1]      [,2]      [,3]
## [1,] 0.2248830 0.4162161 0.5454183
## [2,] 0.2903342 0.5791370 0.2187618

2.7 Tensor Product Operations

There are some tensor specific product such as Hadamard product, Kronecker product, and Khatri-Rao product.

2.7.1 Hadamard 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\).

Figure 7: Hadamard Product

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

2.7.2 Kronecker Product

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}\).

Figure 8: Kronecker Product

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

2.7.3 Khatri-Rao Product

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} \]

Figure 9: Khatri-Rao Product

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

2.8 Utilities Functions

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] -852642276 1261929587
##   .. .. .. ..@ 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] -852642276 1261929587
##   .. .. .. ..@ 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] -852642276 1261929587
##   .. .. .. ..@ sparse  : logi FALSE

2.8.1 Bind Operations

modebind_list collapses multiple DelayedArray objects into single DelayedArray object.

m specifies the collapsed dimension.

Figure 10: Bind Operations

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

Session information

## R version 4.4.0 beta (2024-04-15 r86425)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 22.04.4 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.19-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.0
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_US.UTF-8        LC_COLLATE=en_US.UTF-8    
##  [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       
## 
## time zone: America/New_York
## tzcode source: system (glibc)
## 
## attached base packages:
## [1] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
##  [1] DelayedRandomArray_1.12.0 HDF5Array_1.32.0         
##  [3] rhdf5_2.48.0              DelayedArray_0.30.0      
##  [5] SparseArray_1.4.0         S4Arrays_1.4.0           
##  [7] abind_1.4-5               IRanges_2.38.0           
##  [9] S4Vectors_0.42.0          MatrixGenerics_1.16.0    
## [11] matrixStats_1.3.0         BiocGenerics_0.50.0      
## [13] Matrix_1.7-0              DelayedTensor_1.10.0     
## [15] BiocStyle_2.32.0         
## 
## loaded via a namespace (and not attached):
##  [1] jsonlite_1.8.8      compiler_4.4.0      BiocManager_1.30.22
##  [4] crayon_1.5.2        rsvd_1.0.5          Rcpp_1.0.12        
##  [7] rhdf5filters_1.16.0 parallel_4.4.0      jquerylib_0.1.4    
## [10] BiocParallel_1.38.0 yaml_2.3.8          fastmap_1.1.1      
## [13] lattice_0.22-6      R6_2.5.1            XVector_0.44.0     
## [16] ScaledMatrix_1.12.0 knitr_1.46          einsum_0.1.2       
## [19] bookdown_0.39       bslib_0.7.0         rlang_1.1.3        
## [22] cachem_1.0.8        xfun_0.43           sass_0.4.9         
## [25] cli_3.6.2           Rhdf5lib_1.26.0     BiocSingular_1.20.0
## [28] zlibbioc_1.50.0     digest_0.6.35       grid_4.4.0         
## [31] irlba_2.3.5.1       rTensor_1.4.8       dqrng_0.3.2        
## [34] lifecycle_1.0.4     evaluate_0.23       codetools_0.2-20   
## [37] beachmat_2.20.0     rmarkdown_2.26      tools_4.4.0        
## [40] htmltools_0.5.8.1