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,] 0.7044129 1.0925209 1.1792159   . 0.764944245 1.488692467
## [2,] 0.6748770 0.8215438 0.7354052   . 0.676497854 1.077956836
## [3,] 1.6798836 1.7432628 0.7414468   . 0.467313554 1.267908968
## [4,] 1.1604300 1.2261038 0.3411192   . 0.004120094 0.935269049
## [5,] 0.1182361 0.8632553 0.6876657   . 0.979105625 0.114139350
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,] 0.7044129 1.0925209 1.1792159   . 0.764944245 1.488692467
## [2,] 0.6748770 0.8215438 0.7354052   . 0.676497854 1.077956836
## [3,] 1.6798836 1.7432628 0.7414468   . 0.467313554 1.267908968
## [4,] 1.1604300 1.2261038 0.3411192   . 0.004120094 0.935269049
## [5,] 0.1182361 0.8632553 0.6876657   . 0.979105625 0.114139350

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] 101789
# Compared to a dense array.
tmp2 <- tempfile(fileext=".h5")
out <- HDF5Array::writeHDF5Array(x, tmp2, "data")
file.info(tmp2)[["size"]]
## [1] 279094
# 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.8471489  0.8471489  0.8471489   .  0.8471489  0.8471489
##    [2,]  0.7277868  0.7277868  0.7277868   .  0.7277868  0.7277868
##    [3,]  0.5261113  0.5261113  0.5261113   .  0.5261113  0.5261113
##    [4,]  0.9010345  0.9010345  0.9010345   .  0.9010345  0.9010345
##    [5,]  0.5261555  0.5261555  0.5261555   .  0.5261555  0.5261555
##     ...          .          .          .   .          .          .
##  [996,] 0.31229104 0.31229104 0.31229104   . 0.31229104 0.31229104
##  [997,] 0.80857800 0.80857800 0.80857800   . 0.80857800 0.80857800
##  [998,] 0.16979425 0.16979425 0.16979425   . 0.16979425 0.16979425
##  [999,] 0.72196537 0.72196537 0.72196537   . 0.72196537 0.72196537
## [1000,] 0.04939984 0.04939984 0.04939984   . 0.04939984 0.04939984
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.8471489  0.8471489  0.8471489   .  0.8471489  0.8471489
##    [2,]  0.7277868  0.7277868  0.7277868   .  0.7277868  0.7277868
##    [3,]  0.5261113  0.5261113  0.5261113   .  0.5261113  0.5261113
##    [4,]  0.9010345  0.9010345  0.9010345   .  0.9010345  0.9010345
##    [5,]  0.5261555  0.5261555  0.5261555   .  0.5261555  0.5261555
##     ...          .          .          .   .          .          .
##  [996,] 0.31229104 0.31229104 0.31229104   . 0.31229104 0.31229104
##  [997,] 0.80857800 0.80857800 0.80857800   . 0.80857800 0.80857800
##  [998,] 0.16979425 0.16979425 0.16979425   . 0.16979425 0.16979425
##  [999,] 0.72196537 0.72196537 0.72196537   . 0.72196537 0.72196537
## [1000,] 0.04939984 0.04939984 0.04939984   . 0.04939984 0.04939984

Session information

sessionInfo()
## R version 4.3.0 RC (2023-04-13 r84269 ucrt)
## Platform: x86_64-w64-mingw32/x64 (64-bit)
## Running under: Windows Server 2022 x64 (build 20348)
## 
## Matrix products: default
## 
## 
## locale:
## [1] LC_COLLATE=C                          
## [2] LC_CTYPE=English_United States.utf8   
## [3] LC_MONETARY=English_United States.utf8
## [4] LC_NUMERIC=C                          
## [5] LC_TIME=English_United States.utf8    
## 
## time zone: America/New_York
## tzcode source: internal
## 
## 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