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> DelayedMatrix object of type "double":
##           [,1]      [,2]      [,3] ...       [,9]      [,10]
## [1,] 0.4635056 1.2184827 1.2581445   . 0.25772267 0.43820970
## [2,] 3.2898702 3.1307929 1.9608751   . 1.52644969 3.28565683
## [3,] 0.7689798 1.0841728 0.9894161   . 0.05591575 1.13348024
## [4,] 0.3753422 0.1973420 1.0905360   . 0.95991765 1.08260509
## [5,] 1.1533418 0.6111023 0.2174118   . 1.02121195 0.26774830
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> DelayedMatrix object of type "double":
##           [,1]      [,2]      [,3] ...       [,9]      [,10]
## [1,] 0.4635056 1.2184827 1.2581445   . 0.25772267 0.43820970
## [2,] 3.2898702 3.1307929 1.9608751   . 1.52644969 3.28565683
## [3,] 0.7689798 1.0841728 0.9894161   . 0.05591575 1.13348024
## [4,] 0.3753422 0.1973420 1.0905360   . 0.95991765 1.08260509
## [5,] 1.1533418 0.6111023 0.2174118   . 1.02121195 0.26774830

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] 101967
# Compared to a dense array.
tmp2 <- tempfile(fileext=".h5")
out <- HDF5Array::writeHDF5Array(x, tmp2, "data")
file.info(tmp2)[["size"]]
## [1] 280464
# 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> DelayedMatrix object of type "double":
##               [,1]       [,2]       [,3] ...     [,999]    [,1000]
##    [1,]  0.6648442  0.6648442  0.6648442   .  0.6648442  0.6648442
##    [2,]  0.1601088  0.1601088  0.1601088   .  0.1601088  0.1601088
##    [3,]  0.9921015  0.9921015  0.9921015   .  0.9921015  0.9921015
##    [4,]  0.2470631  0.2470631  0.2470631   .  0.2470631  0.2470631
##    [5,]  0.3865734  0.3865734  0.3865734   .  0.3865734  0.3865734
##     ...          .          .          .   .          .          .
##  [996,] 0.41431143 0.41431143 0.41431143   . 0.41431143 0.41431143
##  [997,] 0.21322712 0.21322712 0.21322712   . 0.21322712 0.21322712
##  [998,] 0.05534847 0.05534847 0.05534847   . 0.05534847 0.05534847
##  [999,] 0.43062090 0.43062090 0.43062090   . 0.43062090 0.43062090
## [1000,] 0.61097723 0.61097723 0.61097723   . 0.61097723 0.61097723
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> DelayedMatrix object of type "double":
##               [,1]       [,2]       [,3] ...     [,999]    [,1000]
##    [1,]  0.6648442  0.6648442  0.6648442   .  0.6648442  0.6648442
##    [2,]  0.1601088  0.1601088  0.1601088   .  0.1601088  0.1601088
##    [3,]  0.9921015  0.9921015  0.9921015   .  0.9921015  0.9921015
##    [4,]  0.2470631  0.2470631  0.2470631   .  0.2470631  0.2470631
##    [5,]  0.3865734  0.3865734  0.3865734   .  0.3865734  0.3865734
##     ...          .          .          .   .          .          .
##  [996,] 0.41431143 0.41431143 0.41431143   . 0.41431143 0.41431143
##  [997,] 0.21322712 0.21322712 0.21322712   . 0.21322712 0.21322712
##  [998,] 0.05534847 0.05534847 0.05534847   . 0.05534847 0.05534847
##  [999,] 0.43062090 0.43062090 0.43062090   . 0.43062090 0.43062090
## [1000,] 0.61097723 0.61097723 0.61097723   . 0.61097723 0.61097723

Session information

sessionInfo()
## R Under development (unstable) (2024-10-21 r87258)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.1 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.21-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.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.35.0      rhdf5_2.51.0          chihaya_1.7.0        
##  [4] DelayedArray_0.33.0   SparseArray_1.7.0     S4Arrays_1.7.0       
##  [7] abind_1.4-8           IRanges_2.41.0        S4Vectors_0.45.0     
## [10] MatrixGenerics_1.19.0 matrixStats_1.4.1     BiocGenerics_0.53.0  
## [13] Matrix_1.7-1          BiocStyle_2.35.0     
## 
## loaded via a namespace (and not attached):
##  [1] crayon_1.5.3        cli_3.6.3           knitr_1.48         
##  [4] rlang_1.1.4         xfun_0.48           jsonlite_1.8.9     
##  [7] htmltools_0.5.8.1   sass_0.4.9          rmarkdown_2.28     
## [10] grid_4.5.0          evaluate_1.0.1      jquerylib_0.1.4    
## [13] fastmap_1.2.0       Rhdf5lib_1.29.0     yaml_2.3.10        
## [16] lifecycle_1.0.4     bookdown_0.41       BiocManager_1.30.25
## [19] compiler_4.5.0      Rcpp_1.0.13         rhdf5filters_1.19.0
## [22] XVector_0.47.0      lattice_0.22-6      digest_0.6.37      
## [25] R6_2.5.1            bslib_0.8.0         tools_4.5.0        
## [28] zlibbioc_1.53.0     cachem_1.1.0