A DelayedArray backend for TileDB

Introduction

TileDB implements a framework for local and remote storage of dense and sparse arrays. We can use this as a DelayedArray backend to provide an array-level abstraction, thus allowing the data to be used in many places where an ordinary array or matrix might be used. The TileDBArray package implements the necessary wrappers around TileDB-R to support read/write operations on TileDB arrays within the DelayedArray framework.

Creating a TileDBArray

Creating a TileDBArray is as easy as:

X <- matrix(rnorm(1000), ncol=10)
library(TileDBArray)
writeTileDBArray(X)
## <100 x 10> TileDBMatrix object of type "double":
##               [,1]        [,2]        [,3] ...        [,9]       [,10]
##   [1,]   1.4132676   0.5368621  -0.6389776   . -0.68125214  0.71990424
##   [2,]  -1.0058837  -0.1903407  -0.2326969   . -0.38328156 -0.24338703
##   [3,]  -2.7218941  -1.6582458   0.3038247   . -0.19806509 -0.83717103
##   [4,]   2.3298569  -1.1237167  -1.5853202   . -0.21934574  0.90496418
##   [5,]   2.0613395  -0.1527372  -0.5673429   .  0.07834146  0.47770955
##    ...           .           .           .   .           .           .
##  [96,] -1.23728377  0.27935716  0.81422601   .  0.10540368 -1.02710371
##  [97,] -0.93764147  1.36130891  0.44619016   .  0.93148089  0.26972302
##  [98,] -0.08979487  2.26682852  0.41272109   . -0.17341739 -0.07249934
##  [99,] -1.63007748 -0.31197688  0.54436964   .  1.17077993  0.22296505
## [100,]  0.56683287  2.02162402  0.18538970   .  0.04430973 -1.87433446

Alternatively, we can use coercion methods:

as(X, "TileDBArray")
## <100 x 10> TileDBMatrix object of type "double":
##               [,1]        [,2]        [,3] ...        [,9]       [,10]
##   [1,]   1.4132676   0.5368621  -0.6389776   . -0.68125214  0.71990424
##   [2,]  -1.0058837  -0.1903407  -0.2326969   . -0.38328156 -0.24338703
##   [3,]  -2.7218941  -1.6582458   0.3038247   . -0.19806509 -0.83717103
##   [4,]   2.3298569  -1.1237167  -1.5853202   . -0.21934574  0.90496418
##   [5,]   2.0613395  -0.1527372  -0.5673429   .  0.07834146  0.47770955
##    ...           .           .           .   .           .           .
##  [96,] -1.23728377  0.27935716  0.81422601   .  0.10540368 -1.02710371
##  [97,] -0.93764147  1.36130891  0.44619016   .  0.93148089  0.26972302
##  [98,] -0.08979487  2.26682852  0.41272109   . -0.17341739 -0.07249934
##  [99,] -1.63007748 -0.31197688  0.54436964   .  1.17077993  0.22296505
## [100,]  0.56683287  2.02162402  0.18538970   .  0.04430973 -1.87433446

This process works also for sparse matrices:

Y <- Matrix::rsparsematrix(1000, 1000, density=0.01)
writeTileDBArray(Y)
## <1000 x 1000> sparse TileDBMatrix object of type "double":
##            [,1]    [,2]    [,3] ...  [,999] [,1000]
##    [1,]    0.00    0.00    0.00   .       0       0
##    [2,]    0.00    0.00    0.00   .       0       0
##    [3,]    0.00    0.00    0.00   .       0       0
##    [4,]    0.00   -0.07    0.00   .       0       0
##    [5,]    0.00    0.00    0.00   .       0       0
##     ...       .       .       .   .       .       .
##  [996,]       0       0       0   .       0       0
##  [997,]       0       0       0   .       0       0
##  [998,]       0       0       0   .       0       0
##  [999,]       0       0       0   .       0       0
## [1000,]       0       0       0   .       0       0

Logical and integer matrices are supported:

writeTileDBArray(Y > 0)
## <1000 x 1000> sparse TileDBMatrix object of type "logical":
##            [,1]    [,2]    [,3] ...  [,999] [,1000]
##    [1,]   FALSE   FALSE   FALSE   .   FALSE   FALSE
##    [2,]   FALSE   FALSE   FALSE   .   FALSE   FALSE
##    [3,]   FALSE   FALSE   FALSE   .   FALSE   FALSE
##    [4,]   FALSE   FALSE   FALSE   .   FALSE   FALSE
##    [5,]   FALSE   FALSE   FALSE   .   FALSE   FALSE
##     ...       .       .       .   .       .       .
##  [996,]   FALSE   FALSE   FALSE   .   FALSE   FALSE
##  [997,]   FALSE   FALSE   FALSE   .   FALSE   FALSE
##  [998,]   FALSE   FALSE   FALSE   .   FALSE   FALSE
##  [999,]   FALSE   FALSE   FALSE   .   FALSE   FALSE
## [1000,]   FALSE   FALSE   FALSE   .   FALSE   FALSE

As are matrices with dimension names:

rownames(X) <- sprintf("GENE_%i", seq_len(nrow(X)))
colnames(X) <- sprintf("SAMP_%i", seq_len(ncol(X)))
writeTileDBArray(X)
## <100 x 10> TileDBMatrix object of type "double":
##               SAMP_1      SAMP_2      SAMP_3 ...      SAMP_9     SAMP_10
##   GENE_1   1.4132676   0.5368621  -0.6389776   . -0.68125214  0.71990424
##   GENE_2  -1.0058837  -0.1903407  -0.2326969   . -0.38328156 -0.24338703
##   GENE_3  -2.7218941  -1.6582458   0.3038247   . -0.19806509 -0.83717103
##   GENE_4   2.3298569  -1.1237167  -1.5853202   . -0.21934574  0.90496418
##   GENE_5   2.0613395  -0.1527372  -0.5673429   .  0.07834146  0.47770955
##      ...           .           .           .   .           .           .
##  GENE_96 -1.23728377  0.27935716  0.81422601   .  0.10540368 -1.02710371
##  GENE_97 -0.93764147  1.36130891  0.44619016   .  0.93148089  0.26972302
##  GENE_98 -0.08979487  2.26682852  0.41272109   . -0.17341739 -0.07249934
##  GENE_99 -1.63007748 -0.31197688  0.54436964   .  1.17077993  0.22296505
## GENE_100  0.56683287  2.02162402  0.18538970   .  0.04430973 -1.87433446

Manipulating TileDBArrays

TileDBArrays are simply DelayedArray objects and can be manipulated as such. The usual conventions for extracting data from matrix-like objects work as expected:

out <- as(X, "TileDBArray")
dim(out)
## [1] 100  10
head(rownames(out))
## [1] "GENE_1" "GENE_2" "GENE_3" "GENE_4" "GENE_5" "GENE_6"
head(out[,1])
##       GENE_1       GENE_2       GENE_3       GENE_4       GENE_5       GENE_6 
##  1.413267588 -1.005883714 -2.721894073  2.329856921  2.061339476 -0.006928295

We can also perform manipulations like subsetting and arithmetic. Note that these operations do not affect the data in the TileDB backend; rather, they are delayed until the values are explicitly required, hence the creation of the DelayedMatrix object.

out[1:5,1:5] 
## <5 x 5> DelayedMatrix object of type "double":
##             SAMP_1      SAMP_2      SAMP_3      SAMP_4      SAMP_5
## GENE_1  1.41326759  0.53686206 -0.63897763 -0.25316857  0.87902877
## GENE_2 -1.00588371 -0.19034067 -0.23269691 -0.71314207 -1.12387684
## GENE_3 -2.72189407 -1.65824584  0.30382470 -0.34977176  0.09407254
## GENE_4  2.32985692 -1.12371671 -1.58532017 -0.22584428 -0.51544422
## GENE_5  2.06133948 -0.15273716 -0.56734294 -1.47620863 -0.14594437
out * 2
## <100 x 10> DelayedMatrix object of type "double":
##              SAMP_1     SAMP_2     SAMP_3 ...      SAMP_9     SAMP_10
##   GENE_1  2.8265352  1.0737241 -1.2779553   .  -1.3625043   1.4398085
##   GENE_2 -2.0117674 -0.3806813 -0.4653938   .  -0.7665631  -0.4867741
##   GENE_3 -5.4437881 -3.3164917  0.6076494   .  -0.3961302  -1.6743421
##   GENE_4  4.6597138 -2.2474334 -3.1706403   .  -0.4386915   1.8099284
##   GENE_5  4.1226790 -0.3054743 -1.1346859   .   0.1566829   0.9554191
##      ...          .          .          .   .           .           .
##  GENE_96 -2.4745675  0.5587143  1.6284520   .  0.21080737 -2.05420743
##  GENE_97 -1.8752829  2.7226178  0.8923803   .  1.86296178  0.53944603
##  GENE_98 -0.1795897  4.5336570  0.8254422   . -0.34683479 -0.14499868
##  GENE_99 -3.2601550 -0.6239538  1.0887393   .  2.34155986  0.44593010
## GENE_100  1.1336657  4.0432480  0.3707794   .  0.08861946 -3.74866893

We can also do more complex matrix operations that are supported by DelayedArray:

colSums(out)
##       SAMP_1       SAMP_2       SAMP_3       SAMP_4       SAMP_5       SAMP_6 
##   3.24065071   0.07319027  10.29918731   1.68261828  -5.54026995 -15.55885952 
##       SAMP_7       SAMP_8       SAMP_9      SAMP_10 
##   4.98271652  -5.53187198  -0.70622341  -4.00485360
out %*% runif(ncol(out))
##                  [,1]
## GENE_1    1.501236265
## GENE_2   -2.269538499
## GENE_3   -0.934427235
## GENE_4   -0.041710955
## GENE_5   -1.290576769
## GENE_6    1.251724230
## GENE_7    2.487364287
## GENE_8   -3.665892649
## GENE_9    1.278634168
## GENE_10   0.344640614
## GENE_11   2.238363679
## GENE_12   0.734710740
## GENE_13  -0.557229756
## GENE_14   0.674444466
## GENE_15  -2.362144364
## GENE_16  -1.655224380
## GENE_17  -0.941097519
## GENE_18   1.077815854
## GENE_19  -0.424557321
## GENE_20  -1.520632481
## GENE_21   1.685243024
## GENE_22  -3.163278762
## GENE_23  -0.871372072
## GENE_24   2.252308160
## GENE_25   2.164595292
## GENE_26  -1.275027083
## GENE_27   0.886615075
## GENE_28   2.240727935
## GENE_29   1.101207975
## GENE_30  -0.685169294
## GENE_31  -2.081056553
## GENE_32  -2.316968592
## GENE_33   1.755509089
## GENE_34  -0.150059230
## GENE_35   0.712533567
## GENE_36   0.998642087
## GENE_37   0.945823864
## GENE_38  -1.176819292
## GENE_39   3.953665317
## GENE_40  -0.625284069
## GENE_41   0.577071967
## GENE_42  -1.198865334
## GENE_43  -0.824307941
## GENE_44  -0.979873238
## GENE_45   1.811838732
## GENE_46  -0.187278764
## GENE_47   0.004256112
## GENE_48  -1.068502748
## GENE_49   4.084897743
## GENE_50   0.513535384
## GENE_51  -3.218591932
## GENE_52  -0.990245435
## GENE_53  -2.088175246
## GENE_54  -1.527269046
## GENE_55  -0.133769346
## GENE_56   0.758515489
## GENE_57   3.400431685
## GENE_58  -0.629575395
## GENE_59   2.410370212
## GENE_60   0.485299355
## GENE_61   2.208002742
## GENE_62   3.192774814
## GENE_63  -0.653397075
## GENE_64  -1.079482422
## GENE_65   1.008422189
## GENE_66  -0.679887469
## GENE_67  -2.657471355
## GENE_68   3.973601698
## GENE_69   3.265722434
## GENE_70  -1.256267256
## GENE_71  -3.181875973
## GENE_72  -0.956983162
## GENE_73  -1.249574665
## GENE_74   0.198575374
## GENE_75  -3.096128181
## GENE_76   1.889497157
## GENE_77  -4.745559032
## GENE_78  -2.369691426
## GENE_79  -0.003484847
## GENE_80   0.010626864
## GENE_81  -2.756907079
## GENE_82   1.732389960
## GENE_83  -1.317125737
## GENE_84   0.377383113
## GENE_85   0.972649215
## GENE_86   1.036523085
## GENE_87  -1.967785128
## GENE_88   0.307906969
## GENE_89  -1.846898774
## GENE_90   0.639740564
## GENE_91   0.375231995
## GENE_92   2.846764374
## GENE_93   0.804657147
## GENE_94   1.213918264
## GENE_95  -2.741227414
## GENE_96  -0.906437928
## GENE_97   1.215887175
## GENE_98  -0.984102269
## GENE_99   1.398566303
## GENE_100  0.254261596

Controlling backend creation

We can adjust some parameters for creating the backend with appropriate arguments to writeTileDBArray(). For example, the example below allows us to control the path to the backend as well as the name of the attribute containing the data.

X <- matrix(rnorm(1000), ncol=10)
path <- tempfile()
writeTileDBArray(X, path=path, attr="WHEE")
## <100 x 10> TileDBMatrix object of type "double":
##               [,1]        [,2]        [,3] ...         [,9]        [,10]
##   [1,]   1.0061774   0.5814132  -1.4886564   .  0.589327662 -0.725769882
##   [2,]   0.7363253   0.4786291   0.1202360   .  0.122558215  0.547715758
##   [3,]   1.0285393  -1.7085957  -0.6931578   . -2.245644533 -0.079730839
##   [4,]  -1.6294840   2.0106551  -0.1657447   .  0.462739609  1.030815958
##   [5,]   0.4092186  -0.5955360  -2.0368308   .  0.002863273  2.645960661
##    ...           .           .           .   .            .            .
##  [96,]  0.48884863 -1.62256125  0.89277515   .    0.5747947   -0.4164272
##  [97,] -0.13991840  1.57946473 -0.22026240   .    0.6547212    0.8963356
##  [98,] -0.20095193  0.42245644 -0.09717781   .   -0.4025164   -1.4839760
##  [99,]  0.05013886 -0.12953153  0.23961957   .   -0.3465574   -0.2498839
## [100,]  1.21002459  0.27104281  1.55494265   .   -0.5296799   -1.1211015

As these arguments cannot be passed during coercion, we instead provide global variables that can be set or unset to affect the outcome.

path2 <- tempfile()
setTileDBPath(path2)
as(X, "TileDBArray") # uses path2 to store the backend.
## <100 x 10> TileDBMatrix object of type "double":
##               [,1]        [,2]        [,3] ...         [,9]        [,10]
##   [1,]   1.0061774   0.5814132  -1.4886564   .  0.589327662 -0.725769882
##   [2,]   0.7363253   0.4786291   0.1202360   .  0.122558215  0.547715758
##   [3,]   1.0285393  -1.7085957  -0.6931578   . -2.245644533 -0.079730839
##   [4,]  -1.6294840   2.0106551  -0.1657447   .  0.462739609  1.030815958
##   [5,]   0.4092186  -0.5955360  -2.0368308   .  0.002863273  2.645960661
##    ...           .           .           .   .            .            .
##  [96,]  0.48884863 -1.62256125  0.89277515   .    0.5747947   -0.4164272
##  [97,] -0.13991840  1.57946473 -0.22026240   .    0.6547212    0.8963356
##  [98,] -0.20095193  0.42245644 -0.09717781   .   -0.4025164   -1.4839760
##  [99,]  0.05013886 -0.12953153  0.23961957   .   -0.3465574   -0.2498839
## [100,]  1.21002459  0.27104281  1.55494265   .   -0.5296799   -1.1211015

Session information

sessionInfo()
## R version 4.6.0 (2026-04-24)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.4 LTS
## 
## Matrix products: default
## BLAS:   /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3 
## LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so;  LAPACK version 3.12.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: Etc/UTC
## tzcode source: system (glibc)
## 
## attached base packages:
## [1] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
##  [1] RcppSpdlog_0.0.29     TileDBArray_1.22.0    DelayedArray_0.38.2  
##  [4] SparseArray_1.12.2    S4Arrays_1.12.0       IRanges_2.46.0       
##  [7] abind_1.4-8           S4Vectors_0.50.1      MatrixGenerics_1.24.0
## [10] matrixStats_1.5.0     BiocGenerics_0.58.1   generics_0.1.4       
## [13] Matrix_1.7-5          BiocStyle_2.40.0     
## 
## loaded via a namespace (and not attached):
##  [1] bit_4.6.0           jsonlite_2.0.0      compiler_4.6.0     
##  [4] BiocManager_1.30.27 Rcpp_1.1.1-1.1      nanoarrow_0.8.0    
##  [7] jquerylib_0.1.4     yaml_2.3.12         fastmap_1.2.0      
## [10] lattice_0.22-9      RcppCCTZ_0.2.14     R6_2.6.1           
## [13] XVector_0.52.0      tiledb_0.33.0       knitr_1.51         
## [16] maketools_1.3.2     bslib_0.11.0        rlang_1.2.0        
## [19] cachem_1.1.0        xfun_0.57           sass_0.4.10        
## [22] sys_3.4.3           bit64_4.8.2         cli_3.6.6          
## [25] spdl_0.0.5          digest_0.6.39       grid_4.6.0         
## [28] lifecycle_1.0.5     data.table_1.18.4   evaluate_1.0.5     
## [31] nanotime_0.3.15     zoo_1.8-15          buildtools_1.0.0   
## [34] rmarkdown_2.31      tools_4.6.0         htmltools_0.5.9