Contents

1 Motivation

The chihaya package saves DelayedArray objects for efficient, portable and stable reproduction of delayed operations in a new R session or other programming frameworks.

Check out the specification for more details.

2 Quick start

Make a DelayedArray object with some operations:

library(DelayedArray)
x <- DelayedArray(matrix(runif(1000), ncol=10))
x <- x[11:15,] / runif(5) 
x <- log2(x + 1)
x
## <5 x 10> matrix of class DelayedMatrix and type "double":
##           [,1]      [,2]      [,3] ...      [,9]     [,10]
## [1,] 1.3055207 0.6368694 1.2564862   . 0.4933343 1.0645851
## [2,] 1.3737801 1.0511653 0.7284169   . 1.0967000 0.6118097
## [3,] 0.6863176 0.8787785 0.7790733   . 0.3551363 0.6856166
## [4,] 0.1659832 1.3458866 1.6040235   . 0.9000644 1.2761396
## [5,] 1.6767594 0.6061622 1.7565022   . 0.8673755 0.8133542
showtree(x)
## 5x10 double: DelayedMatrix object
## └─ 5x10 double: Stack of 2 unary iso op(s)
##    └─ 5x10 double: Unary iso op with args
##       └─ 5x10 double: Subset
##          └─ 100x10 double: [seed] matrix object

Save it into a HDF5 file with saveDelayed():

library(chihaya)
tmp <- tempfile(fileext=".h5")
saveDelayed(x, tmp)
rhdf5::h5ls(tmp)
##                            group    name       otype  dclass      dim
## 0                              / delayed   H5I_GROUP                 
## 1                       /delayed    base H5I_DATASET   FLOAT    ( 0 )
## 2                       /delayed  method H5I_DATASET  STRING    ( 0 )
## 3                       /delayed    seed   H5I_GROUP                 
## 4                  /delayed/seed  method H5I_DATASET  STRING    ( 0 )
## 5                  /delayed/seed    seed   H5I_GROUP                 
## 6             /delayed/seed/seed   along H5I_DATASET INTEGER    ( 0 )
## 7             /delayed/seed/seed  method H5I_DATASET  STRING    ( 0 )
## 8             /delayed/seed/seed    seed   H5I_GROUP                 
## 9        /delayed/seed/seed/seed   index   H5I_GROUP                 
## 10 /delayed/seed/seed/seed/index       0 H5I_DATASET INTEGER        5
## 11       /delayed/seed/seed/seed    seed   H5I_GROUP                 
## 12  /delayed/seed/seed/seed/seed    data H5I_DATASET   FLOAT 100 x 10
## 13  /delayed/seed/seed/seed/seed  native H5I_DATASET INTEGER    ( 0 )
## 14            /delayed/seed/seed    side H5I_DATASET  STRING    ( 0 )
## 15            /delayed/seed/seed   value H5I_DATASET   FLOAT        5
## 16                 /delayed/seed    side H5I_DATASET  STRING    ( 0 )
## 17                 /delayed/seed   value H5I_DATASET   FLOAT    ( 0 )

And then load it back in later:

y <- loadDelayed(tmp)
y
## <5 x 10> matrix of class DelayedMatrix and type "double":
##           [,1]      [,2]      [,3] ...      [,9]     [,10]
## [1,] 1.3055207 0.6368694 1.2564862   . 0.4933343 1.0645851
## [2,] 1.3737801 1.0511653 0.7284169   . 1.0967000 0.6118097
## [3,] 0.6863176 0.8787785 0.7790733   . 0.3551363 0.6856166
## [4,] 0.1659832 1.3458866 1.6040235   . 0.9000644 1.2761396
## [5,] 1.6767594 0.6061622 1.7565022   . 0.8673755 0.8133542

Of course, this is not a particularly interesting case as we end up saving the original array inside our HDF5 file anyway. The real fun begins when you have some more interesting seeds.

3 More interesting seeds

We can use the delayed nature of the operations to avoid breaking sparsity. For example:

library(Matrix)
x <- rsparsematrix(1000, 1000, density=0.01)
x <- DelayedArray(x) + runif(1000)

tmp <- tempfile(fileext=".h5")
saveDelayed(x, tmp)
rhdf5::h5ls(tmp)
##            group     name       otype  dclass   dim
## 0              /  delayed   H5I_GROUP              
## 1       /delayed    along H5I_DATASET INTEGER ( 0 )
## 2       /delayed   method H5I_DATASET  STRING ( 0 )
## 3       /delayed     seed   H5I_GROUP              
## 4  /delayed/seed     data H5I_DATASET   FLOAT 10000
## 5  /delayed/seed dimnames   H5I_GROUP              
## 6  /delayed/seed  indices H5I_DATASET INTEGER 10000
## 7  /delayed/seed   indptr H5I_DATASET INTEGER  1001
## 8  /delayed/seed    shape H5I_DATASET INTEGER     2
## 9       /delayed     side H5I_DATASET  STRING ( 0 )
## 10      /delayed    value H5I_DATASET   FLOAT  1000
file.info(tmp)[["size"]]
## [1] 102114
# Compared to a dense array.
tmp2 <- tempfile(fileext=".h5")
out <- HDF5Array::writeHDF5Array(x, tmp2, "data")
file.info(tmp2)[["size"]]
## [1] 279473
# Loading it back in.
y <- loadDelayed(tmp)
showtree(y)
## 1000x1000 double: DelayedMatrix object
## └─ 1000x1000 double: Unary iso op with args
##    └─ 1000x1000 double, sparse: [seed] dgCMatrix object

We can also store references to external files, thus avoiding data duplication:

library(HDF5Array)
test <- HDF5Array(tmp2, "data")
stuff <- log2(test + 1)
stuff
## <1000 x 1000> matrix of class DelayedMatrix and type "double":
##              [,1]      [,2]      [,3] ...    [,999]   [,1000]
##    [1,] 0.8456612 0.8456612 1.1680567   . 1.0752507 0.8456612
##    [2,] 0.3623056 0.3623056 0.3623056   . 0.3623056 0.3623056
##    [3,] 0.7545169 0.7545169 0.7545169   . 0.7545169 0.7545169
##    [4,] 0.8953220 0.8953220 0.8953220   . 0.8953220 0.8953220
##    [5,] 0.3593807 0.3593807 0.3593807   . 0.3593807 0.3593807
##     ...         .         .         .   .         .         .
##  [996,] 0.6921004 0.6921004 0.6921004   . 0.6921004 0.6921004
##  [997,] 0.8573638 0.8573638 0.8573638   . 0.8573638 0.8573638
##  [998,] 0.5601155 0.5601155 0.5601155   . 0.5601155 0.5601155
##  [999,] 0.8231144 0.8231144 0.8231144   . 0.8231144 0.8231144
## [1000,] 0.7176613 0.7176613 0.7176613   . 0.7176613 0.7176613
tmp <- tempfile(fileext=".h5")
saveDelayed(stuff, tmp)
rhdf5::h5ls(tmp)
##                 group       name       otype  dclass   dim
## 0                   /    delayed   H5I_GROUP              
## 1            /delayed       base H5I_DATASET   FLOAT ( 0 )
## 2            /delayed     method H5I_DATASET  STRING ( 0 )
## 3            /delayed       seed   H5I_GROUP              
## 4       /delayed/seed     method H5I_DATASET  STRING ( 0 )
## 5       /delayed/seed       seed   H5I_GROUP              
## 6  /delayed/seed/seed dimensions H5I_DATASET INTEGER     2
## 7  /delayed/seed/seed       file H5I_DATASET  STRING ( 0 )
## 8  /delayed/seed/seed       name H5I_DATASET  STRING ( 0 )
## 9  /delayed/seed/seed     sparse H5I_DATASET INTEGER ( 0 )
## 10 /delayed/seed/seed       type H5I_DATASET  STRING ( 0 )
## 11      /delayed/seed       side H5I_DATASET  STRING ( 0 )
## 12      /delayed/seed      value H5I_DATASET   FLOAT ( 0 )
file.info(tmp)[["size"]] # size of the delayed operations + pointer to the actual file
## [1] 49642
y <- loadDelayed(tmp)
y
## <1000 x 1000> matrix of class DelayedMatrix and type "double":
##              [,1]      [,2]      [,3] ...    [,999]   [,1000]
##    [1,] 0.8456612 0.8456612 1.1680567   . 1.0752507 0.8456612
##    [2,] 0.3623056 0.3623056 0.3623056   . 0.3623056 0.3623056
##    [3,] 0.7545169 0.7545169 0.7545169   . 0.7545169 0.7545169
##    [4,] 0.8953220 0.8953220 0.8953220   . 0.8953220 0.8953220
##    [5,] 0.3593807 0.3593807 0.3593807   . 0.3593807 0.3593807
##     ...         .         .         .   .         .         .
##  [996,] 0.6921004 0.6921004 0.6921004   . 0.6921004 0.6921004
##  [997,] 0.8573638 0.8573638 0.8573638   . 0.8573638 0.8573638
##  [998,] 0.5601155 0.5601155 0.5601155   . 0.5601155 0.5601155
##  [999,] 0.8231144 0.8231144 0.8231144   . 0.8231144 0.8231144
## [1000,] 0.7176613 0.7176613 0.7176613   . 0.7176613 0.7176613

Session information

sessionInfo()
## R version 4.3.0 RC (2023-04-13 r84269)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 22.04.2 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.17-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_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       
## 
## 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] HDF5Array_1.28.0      rhdf5_2.44.0          chihaya_1.0.0        
##  [4] DelayedArray_0.26.0   IRanges_2.34.0        S4Vectors_0.38.0     
##  [7] MatrixGenerics_1.12.0 matrixStats_0.63.0    BiocGenerics_0.46.0  
## [10] Matrix_1.5-4          BiocStyle_2.28.0     
## 
## loaded via a namespace (and not attached):
##  [1] cli_3.6.1           knitr_1.42          rlang_1.1.0        
##  [4] xfun_0.39           jsonlite_1.8.4      htmltools_0.5.5    
##  [7] sass_0.4.5          rmarkdown_2.21      grid_4.3.0         
## [10] evaluate_0.20       jquerylib_0.1.4     fastmap_1.1.1      
## [13] Rhdf5lib_1.22.0     yaml_2.3.7          bookdown_0.33      
## [16] BiocManager_1.30.20 compiler_4.3.0      Rcpp_1.0.10        
## [19] rhdf5filters_1.12.0 lattice_0.21-8      digest_0.6.31      
## [22] R6_2.5.1            bslib_0.4.2         tools_4.3.0        
## [25] cachem_1.0.7