Authors: Koki Tsuyuzaki [aut, cre]
Last modified: 2022-04-26 14:34:11
Compiled: Tue Apr 26 16:27:03 2022

1 Setting

suppressPackageStartupMessages(library("DelayedTensor"))
suppressPackageStartupMessages(library("DelayedArray"))
suppressPackageStartupMessages(library("HDF5Array"))
suppressPackageStartupMessages(library("DelayedRandomArray"))

darr <- RandomUnifArray(c(3,4,5))

setVerbose(FALSE)
setSparse(FALSE)
setAutoBlockSize(1E+8)
## automatic block size set to 1e+08 bytes (was 1e+08)
tmpdir <- tempdir()
setHDF5DumpDir(tmpdir)

2 Tensor Decomposition

Tensor decomposition models decompose multiple factor matrices and core tensor. Each factor matrix means the patterns of each mode and is used for the visualization and the downstream analysis. Core tensor means the intensity of the patterns and is used to decide which patterns are informative.

We reimplemented some of the tensor decomposition functions of rTensor using block processing of DelayedArray.

Only tensor decomposition algorithms and utility functions that require Fast Fourier Transform (e.g., t_mult, t_svd, and t_svd_reconstruct) are exceptions and have not yet been implemented in DelayedArray because we are still investigating how to calculate them with out-of-core manner.

2.1 Tucker Decomposition

Suppose a tensor \(\mathcal{X} \in \Re^{I \times J \times K}\). Tucker decomposition models decomposes a tensor \(\mathcal{X}\) into a core tensor \(\mathcal{G} \in \Re^{p \times q \times r}\), and multiple factor matrices \(A \in \Re^{I \times p}\), \(B \in \Re^{J \times q}\), and \(C \in \Re^{K \times r}\) (\(p \leq I, q \leq J, \textrm{and}\ r \leq K\)).

\[ \mathcal{X} = \mathcal{G} \times_{1} A \times_{2} B \times_{3} C \]

For simplicity, here we will use a third-order tensor but the Tucker decomposition can be applied to tensors of larger orders.

There are well-known two algorithms; Higher-Order Singular Value Decomposition (HOSVD) and Higher-order Orthogonal Iteration (HOOI).

hosvd performs the HOSVD Tucker decomposition and tucker performs HOOI Tucker decomposition.

For the details, check the hosvd and tucker functions of rTensor.

out_hosvd <- hosvd(darr, ranks=c(2,1,3))
## 
  |                                                                            
  |                                                                      |   0%
  |                                                                            
  |=======================                                               |  33%
  |                                                                            
  |===============================================                       |  67%
  |                                                                            
  |======================================================================| 100%
str(out_hosvd)
## List of 4
##  $ Z          :Formal class 'DelayedArray' [package "DelayedArray"] with 1 slot
##   .. ..@ seed:Formal class 'HDF5ArraySeed' [package "HDF5Array"] with 7 slots
##   .. .. .. ..@ filepath : chr "/tmp/RtmpreF1Pi/auto00125.h5"
##   .. .. .. ..@ name     : chr "/HDF5ArrayAUTO00125"
##   .. .. .. ..@ as_sparse: logi FALSE
##   .. .. .. ..@ type     : chr NA
##   .. .. .. ..@ dim      : int [1:3] 2 1 3
##   .. .. .. ..@ chunkdim : int [1:3] 2 1 3
##   .. .. .. ..@ first_val: num -4.39
##  $ U          :List of 3
##   ..$ :Formal class 'DelayedMatrix' [package "DelayedArray"] with 1 slot
##   .. .. ..@ seed: num [1:3, 1:2] -0.639 -0.581 -0.504 0.384 0.327 ...
##   ..$ :Formal class 'DelayedMatrix' [package "DelayedArray"] with 1 slot
##   .. .. ..@ seed: num [1:4, 1] -0.483 -0.396 -0.575 -0.528
##   ..$ :Formal class 'DelayedMatrix' [package "DelayedArray"] with 1 slot
##   .. .. ..@ seed: num [1:5, 1:3] -0.391 -0.569 -0.417 -0.458 -0.375 ...
##  $ est        :Formal class 'DelayedArray' [package "DelayedArray"] with 1 slot
##   .. ..@ seed:Formal class 'HDF5ArraySeed' [package "HDF5Array"] with 7 slots
##   .. .. .. ..@ filepath : chr "/tmp/RtmpreF1Pi/auto00142.h5"
##   .. .. .. ..@ name     : chr "/HDF5ArrayAUTO00142"
##   .. .. .. ..@ as_sparse: logi FALSE
##   .. .. .. ..@ type     : chr NA
##   .. .. .. ..@ dim      : int [1:3] 3 4 5
##   .. .. .. ..@ chunkdim : int [1:3] 3 4 5
##   .. .. .. ..@ first_val: num 0.51
##  $ fnorm_resid: num 2.05
out_tucker <- tucker(darr, ranks=c(2,3,2))
## 
  |                                                                            
  |                                                                      |   0%
  |                                                                            
  |===                                                                   |   4%
  |                                                                            
  |======                                                                |   8%
  |                                                                            
  |========                                                              |  12%
  |                                                                            
  |===========                                                           |  16%
  |                                                                            
  |==============                                                        |  20%
  |                                                                            
  |=================                                                     |  24%
  |                                                                            
  |====================                                                  |  28%
  |                                                                            
  |======================                                                |  32%
  |                                                                            
  |=========================                                             |  36%
  |                                                                            
  |============================                                          |  40%
  |                                                                            
  |======================================================================| 100%
str(out_tucker)
## List of 7
##  $ Z           :Formal class 'DelayedArray' [package "DelayedArray"] with 1 slot
##   .. ..@ seed:Formal class 'HDF5ArraySeed' [package "HDF5Array"] with 7 slots
##   .. .. .. ..@ filepath : chr "/tmp/RtmpreF1Pi/auto00761.h5"
##   .. .. .. ..@ name     : chr "/HDF5ArrayAUTO00761"
##   .. .. .. ..@ as_sparse: logi FALSE
##   .. .. .. ..@ type     : chr NA
##   .. .. .. ..@ dim      : int [1:3] 2 3 2
##   .. .. .. ..@ chunkdim : int [1:3] 2 3 2
##   .. .. .. ..@ first_val: num -4.39
##  $ U           :List of 3
##   ..$ :Formal class 'DelayedMatrix' [package "DelayedArray"] with 1 slot
##   .. .. ..@ seed: num [1:3, 1:2] -0.633 -0.587 -0.504 0.136 0.557 ...
##   ..$ :Formal class 'DelayedMatrix' [package "DelayedArray"] with 1 slot
##   .. .. ..@ seed: num [1:4, 1:3] -0.481 -0.403 -0.572 -0.528 -0.195 ...
##   ..$ :Formal class 'DelayedMatrix' [package "DelayedArray"] with 1 slot
##   .. .. ..@ seed: num [1:5, 1:2] -0.392 -0.57 -0.402 -0.458 -0.387 ...
##  $ conv        : logi TRUE
##  $ est         :Formal class 'DelayedArray' [package "DelayedArray"] with 1 slot
##   .. ..@ seed:Formal class 'HDF5ArraySeed' [package "HDF5Array"] with 7 slots
##   .. .. .. ..@ filepath : chr "/tmp/RtmpreF1Pi/auto00795.h5"
##   .. .. .. ..@ name     : chr "/HDF5ArrayAUTO00795"
##   .. .. .. ..@ as_sparse: logi FALSE
##   .. .. .. ..@ type     : chr NA
##   .. .. .. ..@ dim      : int [1:3] 3 4 5
##   .. .. .. ..@ chunkdim : int [1:3] 3 4 5
##   .. .. .. ..@ first_val: num 0.409
##  $ norm_percent: num 66
##  $ fnorm_resid : num 1.66
##  $ all_resids  : num [1:10] 1.66 1.66 1.66 1.66 1.66 ...

2.2 CANDECOMP/PARAFAC (CP) Decomposition

Suppose a tensor \(\mathcal{X} \in \Re^{I \times J \times K}\). CP decomposition models decomposes a tensor \(X\) into a core diagonal tensor \(\mathcal{G} \in \Re^{r \times r \times r}\), and multiple factor matrices \(A \in \Re^{I \times r}\), \(B \in \Re^{J \times r}\), and \(C \in \Re^{K \times r}\) (\(r \leq \min{\left(I,J,K\right)}\)).

\[ \mathcal{X} = \mathcal{G} \times_{1} A \times_{2} B \times_{3} C \]

For simplicity, here we will use a third-order tensor but the CP decomposition can be applied to tensors of larger orders.

Alternating least squares (ALS) is a well-known algorithm for CP decomposition and cp performs the ALS CP decomposition.

For the details, check the cp function of rTensor.

out_cp <- cp(darr, num_components=2)
## 
  |                                                                            
  |                                                                      |   0%
  |                                                                            
  |===                                                                   |   4%
  |                                                                            
  |======                                                                |   8%
  |                                                                            
  |========                                                              |  12%
  |                                                                            
  |===========                                                           |  16%
  |                                                                            
  |==============                                                        |  20%
  |                                                                            
  |=================                                                     |  24%
  |                                                                            
  |====================                                                  |  28%
  |                                                                            
  |======================                                                |  32%
  |                                                                            
  |=========================                                             |  36%
  |                                                                            
  |============================                                          |  40%
  |                                                                            
  |===============================                                       |  44%
  |                                                                            
  |==================================                                    |  48%
  |                                                                            
  |====================================                                  |  52%
  |                                                                            
  |=======================================                               |  56%
  |                                                                            
  |==========================================                            |  60%
  |                                                                            
  |=============================================                         |  64%
  |                                                                            
  |================================================                      |  68%
  |                                                                            
  |==================================================                    |  72%
  |                                                                            
  |=====================================================                 |  76%
  |                                                                            
  |========================================================              |  80%
  |                                                                            
  |===========================================================           |  84%
  |                                                                            
  |==============================================================        |  88%
  |                                                                            
  |================================================================      |  92%
  |                                                                            
  |===================================================================   |  96%
  |                                                                            
  |======================================================================| 100%
str(out_cp)
## List of 7
##  $ lambdas     : num [1:2] 8.51 26.1
##  $ U           :List of 3
##   ..$ :Formal class 'DelayedMatrix' [package "DelayedArray"] with 1 slot
##   .. .. ..@ seed:Formal class 'DelayedAperm' [package "DelayedArray"] with 2 slots
##   .. .. .. .. ..@ perm: int [1:2] 2 1
##   .. .. .. .. ..@ seed:Formal class 'DelayedUnaryIsoOpWithArgs' [package "DelayedArray"] with 6 slots
##   .. .. .. .. .. .. ..@ OP    :function (e1, e2)  
##   .. .. .. .. .. .. ..@ Largs : list()
##   .. .. .. .. .. .. ..@ Rargs :List of 1
##   .. .. .. .. .. .. .. ..$ : num [1:2] 8.51 26.12
##   .. .. .. .. .. .. ..@ Lalong: int(0) 
##   .. .. .. .. .. .. ..@ Ralong: int 1
##   .. .. .. .. .. .. ..@ seed  :Formal class 'DelayedAperm' [package "DelayedArray"] with 2 slots
##   .. .. .. .. .. .. .. .. ..@ perm: int [1:2] 2 1
##   .. .. .. .. .. .. .. .. ..@ seed:Formal class 'HDF5ArraySeed' [package "HDF5Array"] with 7 slots
##   .. .. .. .. .. .. .. .. .. .. ..@ filepath : chr "/tmp/RtmpreF1Pi/auto02394.h5"
##   .. .. .. .. .. .. .. .. .. .. ..@ name     : chr "/HDF5ArrayAUTO02394"
##   .. .. .. .. .. .. .. .. .. .. ..@ as_sparse: logi FALSE
##   .. .. .. .. .. .. .. .. .. .. ..@ type     : chr NA
##   .. .. .. .. .. .. .. .. .. .. ..@ dim      : int [1:2] 3 2
##   .. .. .. .. .. .. .. .. .. .. ..@ chunkdim : int [1:2] 3 2
##   .. .. .. .. .. .. .. .. .. .. ..@ first_val: num 0.0206
##   ..$ :Formal class 'DelayedMatrix' [package "DelayedArray"] with 1 slot
##   .. .. ..@ seed:Formal class 'DelayedAperm' [package "DelayedArray"] with 2 slots
##   .. .. .. .. ..@ perm: int [1:2] 2 1
##   .. .. .. .. ..@ seed:Formal class 'DelayedUnaryIsoOpWithArgs' [package "DelayedArray"] with 6 slots
##   .. .. .. .. .. .. ..@ OP    :function (e1, e2)  
##   .. .. .. .. .. .. ..@ Largs : list()
##   .. .. .. .. .. .. ..@ Rargs :List of 1
##   .. .. .. .. .. .. .. ..$ : num [1:2] 8.51 26.1
##   .. .. .. .. .. .. ..@ Lalong: int(0) 
##   .. .. .. .. .. .. ..@ Ralong: int 1
##   .. .. .. .. .. .. ..@ seed  :Formal class 'DelayedAperm' [package "DelayedArray"] with 2 slots
##   .. .. .. .. .. .. .. .. ..@ perm: int [1:2] 2 1
##   .. .. .. .. .. .. .. .. ..@ seed:Formal class 'HDF5ArraySeed' [package "HDF5Array"] with 7 slots
##   .. .. .. .. .. .. .. .. .. .. ..@ filepath : chr "/tmp/RtmpreF1Pi/auto02417.h5"
##   .. .. .. .. .. .. .. .. .. .. ..@ name     : chr "/HDF5ArrayAUTO02417"
##   .. .. .. .. .. .. .. .. .. .. ..@ as_sparse: logi FALSE
##   .. .. .. .. .. .. .. .. .. .. ..@ type     : chr NA
##   .. .. .. .. .. .. .. .. .. .. ..@ dim      : int [1:2] 4 2
##   .. .. .. .. .. .. .. .. .. .. ..@ chunkdim : int [1:2] 4 2
##   .. .. .. .. .. .. .. .. .. .. ..@ first_val: num 3.01
##   ..$ :Formal class 'DelayedMatrix' [package "DelayedArray"] with 1 slot
##   .. .. ..@ seed:Formal class 'DelayedAperm' [package "DelayedArray"] with 2 slots
##   .. .. .. .. ..@ perm: int [1:2] 2 1
##   .. .. .. .. ..@ seed:Formal class 'DelayedUnaryIsoOpWithArgs' [package "DelayedArray"] with 6 slots
##   .. .. .. .. .. .. ..@ OP    :function (e1, e2)  
##   .. .. .. .. .. .. ..@ Largs : list()
##   .. .. .. .. .. .. ..@ Rargs :List of 1
##   .. .. .. .. .. .. .. ..$ : num [1:2] 8.51 26.1
##   .. .. .. .. .. .. ..@ Lalong: int(0) 
##   .. .. .. .. .. .. ..@ Ralong: int 1
##   .. .. .. .. .. .. ..@ seed  :Formal class 'DelayedAperm' [package "DelayedArray"] with 2 slots
##   .. .. .. .. .. .. .. .. ..@ perm: int [1:2] 2 1
##   .. .. .. .. .. .. .. .. ..@ seed:Formal class 'HDF5ArraySeed' [package "HDF5Array"] with 7 slots
##   .. .. .. .. .. .. .. .. .. .. ..@ filepath : chr "/tmp/RtmpreF1Pi/auto02440.h5"
##   .. .. .. .. .. .. .. .. .. .. ..@ name     : chr "/HDF5ArrayAUTO02440"
##   .. .. .. .. .. .. .. .. .. .. ..@ as_sparse: logi FALSE
##   .. .. .. .. .. .. .. .. .. .. ..@ type     : chr NA
##   .. .. .. .. .. .. .. .. .. .. ..@ dim      : int [1:2] 5 2
##   .. .. .. .. .. .. .. .. .. .. ..@ chunkdim : int [1:2] 5 2
##   .. .. .. .. .. .. .. .. .. .. ..@ first_val: num -0.775
##  $ conv        : logi FALSE
##  $ est         :Formal class 'DelayedArray' [package "DelayedArray"] with 1 slot
##   .. ..@ seed:Formal class 'HDF5ArraySeed' [package "HDF5Array"] with 7 slots
##   .. .. .. ..@ filepath : chr "/tmp/RtmpreF1Pi/auto02457.h5"
##   .. .. .. ..@ name     : chr "/HDF5ArrayAUTO02457"
##   .. .. .. ..@ as_sparse: logi FALSE
##   .. .. .. ..@ type     : chr NA
##   .. .. .. ..@ dim      : int [1:3] 3 4 5
##   .. .. .. ..@ chunkdim : int [1:3] 3 4 5
##   .. .. .. ..@ first_val: num 0.466
##  $ norm_percent: num 63.3
##  $ fnorm_resid : num 1.79
##  $ all_resids  : num [1:24] 2.14 2.02 1.86 1.8 1.79 ...

2.3 Multilinear Principal Component Analysis (MPCA)

MPCA is a kind of Tucker decomposition and when the order of tensor is 3, this is also known as the Generalized Low-Rank Approximation of Matrices (GLRAM). For the details, check the mpca function of rTensor.

out_mpca <- mpca(darr, ranks=c(2,2))
## 
  |                                                                            
  |                                                                      |   0%
  |                                                                            
  |===                                                                   |   4%
  |                                                                            
  |======                                                                |   8%
  |                                                                            
  |========                                                              |  12%
  |                                                                            
  |===========                                                           |  16%
  |                                                                            
  |==============                                                        |  20%
  |                                                                            
  |=================                                                     |  24%
  |                                                                            
  |======================================================================| 100%
str(out_mpca)
## List of 7
##  $ Z_ext       :Formal class 'DelayedArray' [package "DelayedArray"] with 1 slot
##   .. ..@ seed:Formal class 'HDF5ArraySeed' [package "HDF5Array"] with 7 slots
##   .. .. .. ..@ filepath : chr "/tmp/RtmpreF1Pi/auto02655.h5"
##   .. .. .. ..@ name     : chr "/HDF5ArrayAUTO02655"
##   .. .. .. ..@ as_sparse: logi FALSE
##   .. .. .. ..@ type     : chr NA
##   .. .. .. ..@ dim      : int [1:3] 2 2 5
##   .. .. .. ..@ chunkdim : int [1:3] 2 2 5
##   .. .. .. ..@ first_val: num 1.69
##  $ U           :List of 3
##   ..$ :Formal class 'DelayedMatrix' [package "DelayedArray"] with 1 slot
##   .. .. ..@ seed: num [1:3, 1:2] -0.632 -0.583 -0.51 0.36 -0.804 ...
##   ..$ :Formal class 'DelayedMatrix' [package "DelayedArray"] with 1 slot
##   .. .. ..@ seed: num [1:4, 1:2] -0.473 -0.398 -0.577 -0.534 0.441 ...
##   ..$ : NULL
##  $ conv        : logi TRUE
##  $ est         :Formal class 'DelayedArray' [package "DelayedArray"] with 1 slot
##   .. ..@ seed:Formal class 'HDF5ArraySeed' [package "HDF5Array"] with 7 slots
##   .. .. .. ..@ filepath : chr "/tmp/RtmpreF1Pi/auto02677.h5"
##   .. .. .. ..@ name     : chr "/HDF5ArrayAUTO02677"
##   .. .. .. ..@ as_sparse: logi FALSE
##   .. .. .. ..@ type     : chr NA
##   .. .. .. ..@ dim      : int [1:3] 3 4 5
##   .. .. .. ..@ chunkdim : int [1:3] 3 4 5
##   .. .. .. ..@ first_val: num 0.476
##  $ norm_percent: num 66
##  $ fnorm_resid : num 1.66
##  $ all_resids  : num [1:6] 1.71 1.68 1.66 1.66 1.66 ...

2.4 Population Value Decomposition (PVD)

Suppose a series of 2D data X_{j}, where \(X_{j} \in \Re^{n_{1} \times n_{2}}\), and \(j = 1, \ldots , n_{3}\). PVD models decomposes a tensor \(X_{j}\) into two common factor matrices across all \(X_{j}\), \(P \in \Re^{n_{1} \times r_{1}}\) and \(D \in \Re^{n_{2} \times r_{2}}\), and a \(X_{j}\) specific factor matrix \(V_{j} \in \Re^{r_{1} \times r_{2}}\) (\(r_{1} \leq n_{1}, \textrm{and}\ r_{2} \leq n_{2}\)).

\[ X_{j} = P \times V_{j} \times D \]

For the details, check the pvd function of rTensor.

out_pvd <- pvd(darr,
  uranks=rep(2,dim(darr)[3]), wranks=rep(3,dim(darr)[3]), a=2, b=3)
## 
  |                                                                            
  |                                                                      |   0%
  |                                                                            
  |=========                                                             |  12%
  |                                                                            
  |==================                                                    |  25%
  |                                                                            
  |==========================                                            |  38%
  |                                                                            
  |===================================                                   |  50%
  |                                                                            
  |============================================                          |  62%
  |                                                                            
  |====================================================                  |  75%
  |                                                                            
  |=============================================================         |  88%
  |                                                                            
  |======================================================================| 100%
str(out_pvd)
## List of 6
##  $ P           :Formal class 'DelayedMatrix' [package "DelayedArray"] with 1 slot
##   .. ..@ seed: num [1:3, 1:2] -0.632 -0.551 -0.545 0.528 -0.821 ...
##  $ D           :Formal class 'DelayedMatrix' [package "DelayedArray"] with 1 slot
##   .. ..@ seed:Formal class 'DelayedAperm' [package "DelayedArray"] with 2 slots
##   .. .. .. ..@ perm: int [1:2] 2 1
##   .. .. .. ..@ seed: num [1:4, 1:3] -0.458 -0.401 -0.608 -0.51 0.385 ...
##  $ V           :List of 5
##   ..$ :Formal class 'DelayedMatrix' [package "DelayedArray"] with 1 slot
##   .. .. ..@ seed:Formal class 'HDF5ArraySeed' [package "HDF5Array"] with 7 slots
##   .. .. .. .. ..@ filepath : chr "/tmp/RtmpreF1Pi/auto02711.h5"
##   .. .. .. .. ..@ name     : chr "/HDF5ArrayAUTO02711"
##   .. .. .. .. ..@ as_sparse: logi FALSE
##   .. .. .. .. ..@ type     : chr NA
##   .. .. .. .. ..@ dim      : int [1:2] 2 3
##   .. .. .. .. ..@ chunkdim : int [1:2] 2 3
##   .. .. .. .. ..@ first_val: num 1.7
##   ..$ :Formal class 'DelayedMatrix' [package "DelayedArray"] with 1 slot
##   .. .. ..@ seed:Formal class 'HDF5ArraySeed' [package "HDF5Array"] with 7 slots
##   .. .. .. .. ..@ filepath : chr "/tmp/RtmpreF1Pi/auto02717.h5"
##   .. .. .. .. ..@ name     : chr "/HDF5ArrayAUTO02717"
##   .. .. .. .. ..@ as_sparse: logi FALSE
##   .. .. .. .. ..@ type     : chr NA
##   .. .. .. .. ..@ dim      : int [1:2] 2 3
##   .. .. .. .. ..@ chunkdim : int [1:2] 2 3
##   .. .. .. .. ..@ first_val: num 2.47
##   ..$ :Formal class 'DelayedMatrix' [package "DelayedArray"] with 1 slot
##   .. .. ..@ seed:Formal class 'HDF5ArraySeed' [package "HDF5Array"] with 7 slots
##   .. .. .. .. ..@ filepath : chr "/tmp/RtmpreF1Pi/auto02723.h5"
##   .. .. .. .. ..@ name     : chr "/HDF5ArrayAUTO02723"
##   .. .. .. .. ..@ as_sparse: logi FALSE
##   .. .. .. .. ..@ type     : chr NA
##   .. .. .. .. ..@ dim      : int [1:2] 2 3
##   .. .. .. .. ..@ chunkdim : int [1:2] 2 3
##   .. .. .. .. ..@ first_val: num 1.76
##   ..$ :Formal class 'DelayedMatrix' [package "DelayedArray"] with 1 slot
##   .. .. ..@ seed:Formal class 'HDF5ArraySeed' [package "HDF5Array"] with 7 slots
##   .. .. .. .. ..@ filepath : chr "/tmp/RtmpreF1Pi/auto02729.h5"
##   .. .. .. .. ..@ name     : chr "/HDF5ArrayAUTO02729"
##   .. .. .. .. ..@ as_sparse: logi FALSE
##   .. .. .. .. ..@ type     : chr NA
##   .. .. .. .. ..@ dim      : int [1:2] 2 3
##   .. .. .. .. ..@ chunkdim : int [1:2] 2 3
##   .. .. .. .. ..@ first_val: num 2.02
##   ..$ :Formal class 'DelayedMatrix' [package "DelayedArray"] with 1 slot
##   .. .. ..@ seed:Formal class 'HDF5ArraySeed' [package "HDF5Array"] with 7 slots
##   .. .. .. .. ..@ filepath : chr "/tmp/RtmpreF1Pi/auto02735.h5"
##   .. .. .. .. ..@ name     : chr "/HDF5ArrayAUTO02735"
##   .. .. .. .. ..@ as_sparse: logi FALSE
##   .. .. .. .. ..@ type     : chr NA
##   .. .. .. .. ..@ dim      : int [1:2] 2 3
##   .. .. .. .. ..@ chunkdim : int [1:2] 2 3
##   .. .. .. .. ..@ first_val: num 1.74
##  $ est         :Formal class 'DelayedArray' [package "DelayedArray"] with 1 slot
##   .. ..@ seed: num [1:3, 1:4, 1:5] 0.55 0.303 0.441 0.626 0.266 ...
##  $ norm_percent: num 64.6
##  $ fnorm_resid : num 1.73

Session information

## R version 4.2.0 RC (2022-04-19 r82224)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.4 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.15-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.15-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.4.0 HDF5Array_1.24.0         rhdf5_2.40.0            
##  [4] DelayedArray_0.22.0      IRanges_2.30.0           S4Vectors_0.34.0        
##  [7] MatrixGenerics_1.8.0     matrixStats_0.62.0       BiocGenerics_0.42.0     
## [10] Matrix_1.4-1             DelayedTensor_1.2.0      BiocStyle_2.24.0        
## 
## loaded via a namespace (and not attached):
##  [1] Rcpp_1.0.8.3        rTensor_1.4.8       bslib_0.3.1        
##  [4] compiler_4.2.0      BiocManager_1.30.17 jquerylib_0.1.4    
##  [7] rhdf5filters_1.8.0  tools_4.2.0         digest_0.6.29      
## [10] jsonlite_1.8.0      evaluate_0.15       lattice_0.20-45    
## [13] rlang_1.0.2         cli_3.3.0           parallel_4.2.0     
## [16] yaml_2.3.5          xfun_0.30           fastmap_1.1.0      
## [19] stringr_1.4.0       knitr_1.38          sass_0.4.1         
## [22] grid_4.2.0          R6_2.5.1            BiocParallel_1.30.0
## [25] rmarkdown_2.14      bookdown_0.26       irlba_2.3.5        
## [28] Rhdf5lib_1.18.0     magrittr_2.0.3      BiocSingular_1.12.0
## [31] htmltools_0.5.2     rsvd_1.0.5          beachmat_2.12.0    
## [34] dqrng_0.3.0         ScaledMatrix_1.4.0  stringi_1.7.6      
## [37] einsum_0.1.0