Contents

1 Introduction

SparseArray is an infrastructure package that enables high-performance sparse data representation and manipulation in R. The workhorse of the package is an array-like container that allows efficient in-memory representation of multidimensional sparse data in R.

2 Install and load the package

Use BiocManager::install() to install the SparseArray package:

if (!require("BiocManager", quietly=TRUE))
    install.packages("BiocManager")
BiocManager::install("SparseArray")

Load the package:

library(SparseArray)

3 The SparseArray virtual class and its two concrete subclasses

The package defines the SparseArray virtual class and two concrete subclasses: COO_SparseArray and SVT_SparseArray.

Each subclass uses its own internal representation of the nonzero multidimensional data: the “COO layout” and the “SVT layout”, respectively.

Note that the SparseArray virtual class makes no assumption about the internal representation of the nonzero data, so it could easily be extended by other S4 classes that use a different layout for the nonzero data.

This vignette focuses on the SVT_SparseArray container, which is the most memory-efficient and feature-complete of the two SparseArray subclasses. The COO_SparseArray class is only provided to support some rare use-cases. In other words, using SVT_SparseArray objects is almost always preferred over using COO_SparseArray objects.

4 SVT_SparseArray objects

The SVT_SparseArray container provides an efficient representation of the nonzero multidimensional data via a novel layout called the “SVT layout”.

Note that SVT_SparseArray objects mimic as much as possible the behavior of ordinary matrix or array objects in base R. In particular, they suppport most of the “standard matrix and array API” defined in base R and in the matrixStats package from CRAN.

4.1 Construction

SVT_SparseArray objects can be constructed in many ways. A common way is to start with an empty object and to subassign nonzero values to it:

svt1 <- SVT_SparseArray(dim=c(6, 4))
svt1[c(1:2, 8, 10, 15:17, 24)] <- (1:8)*10L
svt1
## <6 x 4 SparseMatrix> of type "integer" [nzcount=8 (33%)]:
##      [,1] [,2] [,3] [,4]
## [1,]   10    0    0    0
## [2,]   20   30    0    0
## [3,]    0    0   50    0
## [4,]    0   40   60    0
## [5,]    0    0   70    0
## [6,]    0    0    0   80
svt2 <- SVT_SparseArray(dim=5:3)
svt2[c(1:2, 8, 10, 15:17, 20, 24, 40, 56:60)] <- (1:15)*10L
svt2
## <5 x 4 x 3 SparseArray> of type "integer" [nzcount=15 (25%)]:
## ,,1
##      [,1] [,2] [,3] [,4]
## [1,]   10    0    0   60
## [2,]   20    0    0   70
## [3,]    0   30    0    0
## [4,]    0    0    0    0
## [5,]    0   40   50   80
## 
## ,,2
##      [,1] [,2] [,3] [,4]
## [1,]    0    0    0    0
## [2,]    0    0    0    0
## [3,]    0    0    0    0
## [4,]   90    0    0    0
## [5,]    0    0    0  100
## 
## ,,3
##      [,1] [,2] [,3] [,4]
## [1,]    0    0    0  110
## [2,]    0    0    0  120
## [3,]    0    0    0  130
## [4,]    0    0    0  140
## [5,]    0    0    0  150

Another way is to coerce a matrix- or array-like object to SVT_SparseArray:

# Coerce a dgCMatrix object to SVT_SparseArray:
dgcm <- Matrix::rsparsematrix(12, 5, density=0.15)
svt3 <- as(dgcm, "SVT_SparseArray")

# Coerce a TENxMatrix object to SVT_SparseArray:
suppressMessages(library(HDF5Array))
M <- writeTENxMatrix(svt3)
svt3b <- as(M, "SVT_SparseArray")

# Sanity check:
stopifnot(identical(svt3, svt3b))

Alternatively, these coercions can be done by simply passing the object to coerce to the SVT_SparseArray() constructor function:

svt3  <- SVT_SparseArray(dgcm)  # same as as(dgcm, "SVT_SparseArray")
svt3b <- SVT_SparseArray(M)     # same as as(M, "SVT_SparseArray")

See ?SVT_SparseArray for more information about the SVT_SparseArray() constructor function and additional examples.

4.2 SVT_SparseArray vs COO_SparseArray

As mentioned earlier, SVT_SparseArray objects are almost always preferred over using COO_SparseArray objects. Coercing to SparseArray or using the SparseArray() constructor function reflects this preference i.e. in both cases the actual class of the returned SparseArray derivative will almost always be SVT_SparseArray (or SVT_SparseMatrix). Except in the rare situation where returning a COO_SparseArray object is a more natural fit for the input object.

For example coercing the following objects to SparseArray will always produce an SVT_SparseArray object:

# Coerce an ordinary matrix to SparseArray:
a <- array(rpois(80, lambda=0.35), dim=c(5, 8, 2))
class(as(a, "SparseArray"))  # SVT_SparseArray
## [1] "SVT_SparseArray"
## attr(,"package")
## [1] "SparseArray"
# Coerce a dgCMatrix object to SparseArray:
svt3  <- as(dgcm, "SparseArray")
class(svt3)  # SVT_SparseArray
## [1] "SVT_SparseMatrix"
## attr(,"package")
## [1] "SparseArray"
# Coerce a TENxMatrix object to SparseArray:
svt3b <- as(M, "SparseArray")
class(svt3)  # SVT_SparseArray
## [1] "SVT_SparseMatrix"
## attr(,"package")
## [1] "SparseArray"

Also using the SparseArray() constructor function on these objects will always produce an SVT_SparseArray object:

SparseArray(a)              # same as as(a, "SparseArray")
## <5 x 8 x 2 SparseArray> of type "integer" [nzcount=23 (29%)]:
## ,,1
##      [,1] [,2] [,3] [,4] [,5] [,6] [,7] [,8]
## [1,]    0    1    1    0    0    1    1    0
## [2,]    0    0    0    0    1    1    1    0
## [3,]    0    0    0    0    0    0    0    0
## [4,]    1    0    0    1    0    0    0    1
## [5,]    3    0    1    0    0    0    1    0
## 
## ,,2
##      [,1] [,2] [,3] [,4] [,5] [,6] [,7] [,8]
## [1,]    0    1    1    0    1    1    1    0
## [2,]    0    0    0    0    0    0    2    0
## [3,]    0    0    0    0    0    1    0    1
## [4,]    0    0    0    0    1    1    0    0
## [5,]    0    0    0    0    0    0    0    0
svt3  <- SparseArray(dgcm)  # same as as(dgcm, "SparseArray")
svt3b <- SparseArray(M)     # same as as(M, "SparseArray")

This is actually the most convenient way to turn an ordinary matrix or array, or a dgCMatrix object, or a TENxMatrix object, into an SVT_SparseArray object.

One situation where as(x, "SparseArray") or SparseArray(x) will return a COO_SparseArray object is when the input object x is a sparseMatrix derivative that uses a compressed row-oriented representation ("R" representation) instead of the more widely used compressed column-oriented representation ("C" representation):

ngrm <- sparseMatrix(i=c(1, 5, 5, 6), j=c(4, 2, 3, 2), repr="R")
class(ngrm)  # ngRMatrix
## [1] "ngRMatrix"
## attr(,"package")
## [1] "Matrix"
class(SparseArray(ngrm))  # COO_SparseMatrix
## [1] "COO_SparseMatrix"
## attr(,"package")
## [1] "SparseArray"

One way to enforce the SVT_SparseArray representation is to coerce the result to SVT_SparseArray:

svt <- as(SparseArray(ngrm), "SVT_SparseArray")
class(svt)  # SVT_SparseMatrix
## [1] "SVT_SparseMatrix"
## attr(,"package")
## [1] "SparseArray"

Finally, note that coercing back to ordinary matrix or array (dense representation) is supported, although obviously not a good idea if the SparseArray object is big:

as.array(svt1)  # same as as.matrix(svt1)
##      [,1] [,2] [,3] [,4]
## [1,]   10    0    0    0
## [2,]   20   30    0    0
## [3,]    0    0   50    0
## [4,]    0   40   60    0
## [5,]    0    0   70    0
## [6,]    0    0    0   80
as.array(svt2)
## , , 1
## 
##      [,1] [,2] [,3] [,4]
## [1,]   10    0    0   60
## [2,]   20    0    0   70
## [3,]    0   30    0    0
## [4,]    0    0    0    0
## [5,]    0   40   50   80
## 
## , , 2
## 
##      [,1] [,2] [,3] [,4]
## [1,]    0    0    0    0
## [2,]    0    0    0    0
## [3,]    0    0    0    0
## [4,]   90    0    0    0
## [5,]    0    0    0  100
## 
## , , 3
## 
##      [,1] [,2] [,3] [,4]
## [1,]    0    0    0  110
## [2,]    0    0    0  120
## [3,]    0    0    0  130
## [4,]    0    0    0  140
## [5,]    0    0    0  150

5 The SparseArray API

5.1 The core array API

SVT_SparseArray objects support the “core array API” defined in base R:

dim(svt2)
## [1] 5 4 3
length(svt2)
## [1] 60
dimnames(svt2) <- list(NULL, letters[1:4], LETTERS[1:3])
svt2
## <5 x 4 x 3 SparseArray> of type "integer" [nzcount=15 (25%)]:
## ,,A
##       a  b  c  d
## [1,] 10  0  0 60
## [2,] 20  0  0 70
## [3,]  0 30  0  0
## [4,]  0  0  0  0
## [5,]  0 40 50 80
## 
## ,,B
##        a   b   c   d
## [1,]   0   0   0   0
## [2,]   0   0   0   0
## [3,]   0   0   0   0
## [4,]  90   0   0   0
## [5,]   0   0   0 100
## 
## ,,C
##        a   b   c   d
## [1,]   0   0   0 110
## [2,]   0   0   0 120
## [3,]   0   0   0 130
## [4,]   0   0   0 140
## [5,]   0   0   0 150

5.2 type() and is_sparse()

type() and is_sparse() are generic functions defined in BiocGenerics and S4Arrays, respectively. They extend the “core array API” defined in base R:

type(svt1)
## [1] "integer"
type(svt1) <- "double"
svt1
## <6 x 4 SparseMatrix> of type "double" [nzcount=8 (33%)]:
##      [,1] [,2] [,3] [,4]
## [1,]   10    0    0    0
## [2,]   20   30    0    0
## [3,]    0    0   50    0
## [4,]    0   40   60    0
## [5,]    0    0   70    0
## [6,]    0    0    0   80
is_sparse(svt1)
## [1] TRUE

See ?SparseArray for more information and additional examples.

5.3 is_nonzero() and the nz*() functions

A set of functions is provided for direct manipulation of the nonzero array elements:

is_nonzero(svt1)
## <6 x 4 SparseMatrix> of type "logical" [nzcount=8 (33%)]:
##       [,1]  [,2]  [,3]  [,4]
## [1,]  TRUE FALSE FALSE FALSE
## [2,]  TRUE  TRUE FALSE FALSE
## [3,] FALSE FALSE  TRUE FALSE
## [4,] FALSE  TRUE  TRUE FALSE
## [5,] FALSE FALSE  TRUE FALSE
## [6,] FALSE FALSE FALSE  TRUE
## Get the number of nonzero array elements in 'svt1':
nzcount(svt1)
## [1] 8
## Extract the "linear indices" of the nonzero array elements in 'svt1':
nzwhich(svt1)
## [1]  1  2  8 10 15 16 17 24
## Extract the "array indices" (a.k.a. "array coordinates") of the
## nonzero array elements in 'svt1':
nzwhich(svt1, arr.ind=TRUE)
##      [,1] [,2]
## [1,]    1    1
## [2,]    2    1
## [3,]    2    2
## [4,]    4    2
## [5,]    3    3
## [6,]    4    3
## [7,]    5    3
## [8,]    6    4
## Extract the values of the nonzero array elements in 'svt1':
nzvals(svt1)
## [1] 10 20 30 40 50 60 70 80

Note that the vectors produced by nzwhich() and nzvals() are parallel, that is, they have the same length and the i-th element in one vector corresponds to the i-th element in the other vector.

See ?is_nonzero for more information and additional examples.

5.4 Subsetting and subassignment

svt2[5:3, , "C"]
## <3 x 4 SparseMatrix> of type "integer" [nzcount=3 (25%)]:
##        a   b   c   d
## [1,]   0   0   0 150
## [2,]   0   0   0 140
## [3,]   0   0   0 130

Like with ordinary arrays in base R, assigning values of type "double" to an SVT_SparseArray object of type "integer" will automatically change the type of the object to "double":

type(svt2)
## [1] "integer"
svt2[5, 1, 3] <- NaN
type(svt2)
## [1] "double"

See ?SparseArray_subsetting for more information and additional examples.

5.5 Summarization methods (whole array)

The following summarization methods are provided at the moment: anyNA(), any, all, min, max, range, sum, prod, mean, var, sd.

anyNA(svt2)
## [1] TRUE
range(svt2, na.rm=TRUE)
## [1]   0 150
mean(svt2, na.rm=TRUE)
## [1] 20.33898
var(svt2, na.rm=TRUE)
## [1] 1717.124

See ?SparseArray_summarization for more information and additional examples.

5.6 Operations from the Arith, Compare, Logic, Math, Math2, and Complex groups

SVT_SparseArray objects support operations from the Arith, Compare, Logic, Math, Math2, and Complex groups, with some restrictions. See ?S4groupGeneric in the methods package for more information about these group generics.

signif((svt1^1.5 + svt1) %% 100 - 0.6 * svt1, digits=2)
## <6 x 4 SparseMatrix> of type "double" [nzcount=8 (33%)]:
##       [,1]  [,2]  [,3]  [,4]
## [1,]  36.0   0.0   0.0   0.0
## [2,]  -2.6  76.0   0.0   0.0
## [3,]   0.0   0.0 -26.0   0.0
## [4,]   0.0  69.0 -11.0   0.0
## [5,]   0.0   0.0  14.0   0.0
## [6,]   0.0   0.0   0.0  48.0

See ?SparseArray_Arith, ?SparseArray_Compare, ?SparseArray_Logic, ?SparseArray_Math, and ?SparseArray_Complex for more information and additional examples.

6 The 2D API

6.1 SVT_SparseMatrix objects

SVT_SparseMatrix objects are just two-dimensional SVT_SparseArray objects. See ?SparseArray for a diagram of the SparseArray class hierarchy.

6.2 Transposition

t(svt1)
## <4 x 6 SparseMatrix> of type "double" [nzcount=8 (33%)]:
##      [,1] [,2] [,3] [,4] [,5] [,6]
## [1,]   10   20    0    0    0    0
## [2,]    0   30    0   40    0    0
## [3,]    0    0   50   60   70    0
## [4,]    0    0    0    0    0   80

Note that multidimensional transposition is supported via aperm():

aperm(svt2)
## <3 x 4 x 5 SparseArray> of type "double" [nzcount=16 (27%)]:
## ,,1
##     a   b   c   d
## A  10   0   0  60
## B   0   0   0   0
## C   0   0   0 110
## 
## ,,2
##     a   b   c   d
## A  20   0   0  70
## B   0   0   0   0
## C   0   0   0 120
## 
## ,,3
##     a   b   c   d
## A   0  30   0   0
## B   0   0   0   0
## C   0   0   0 130
## 
## ,,4
##     a   b   c   d
## A   0   0   0   0
## B  90   0   0   0
## C   0   0   0 140
## 
## ,,5
##     a   b   c   d
## A   0  40  50  80
## B   0   0   0 100
## C NaN   0   0 150

See ?SparseArray_aperm for more information and additional examples.

6.3 Combine multidimensional objects along a given dimension

Like ordinary matrices in base R, SVT_SparseMatrix objects can be combined by rows or columns, with rbind() or cbind():

svt4 <- poissonSparseMatrix(6, 2, density=0.5)

cbind(svt1, svt4)
## <6 x 6 SparseMatrix> of type "double" [nzcount=16 (44%)]:
##      [,1] [,2] [,3] [,4] [,5] [,6]
## [1,]   10    0    0    0    0    1
## [2,]   20   30    0    0    2    0
## [3,]    0    0   50    0    1    2
## [4,]    0   40   60    0    0    1
## [5,]    0    0   70    0    3    1
## [6,]    0    0    0   80    0    1

Note that multidimensional objects can be combined along any dimension with abind():

svt5a <- poissonSparseArray(c(5, 6, 2), density=0.4)
svt5b <- poissonSparseArray(c(5, 6, 5), density=0.2)
svt5c <- poissonSparseArray(c(5, 6, 4), density=0.2)
abind(svt5a, svt5b, svt5c)
## <5 x 6 x 11 SparseArray> of type "integer" [nzcount=76 (23%)]:
## ,,1
##      [,1] [,2] [,3] [,4] [,5] [,6]
## [1,]    0    1    0    0    0    0
## [2,]    0    0    0    0    1    0
## [3,]    1    1    1    1    0    0
## [4,]    0    3    1    1    0    1
## [5,]    0    1    0    0    1    0
## 
## ...
## 
## ,,11
##      [,1] [,2] [,3] [,4] [,5] [,6]
## [1,]    0    0    0    0    0    0
## [2,]    0    1    0    1    1    0
## [3,]    0    0    0    1    0    0
## [4,]    0    0    0    1    0    0
## [5,]    0    1    0    0    0    0
svt6a <- aperm(svt5a, c(1, 3:2))
svt6b <- aperm(svt5b, c(1, 3:2))
svt6c <- aperm(svt5c, c(1, 3:2))
abind(svt6a, svt6b, svt6c, along=2)
## <5 x 11 x 6 SparseArray> of type "integer" [nzcount=76 (23%)]:
## ,,1
##       [,1]  [,2]  [,3]  [,4] ...  [,8]  [,9] [,10] [,11]
## [1,]     0     0     0     0   .     0     0     0     0
## [2,]     0     0     0     1   .     0     1     0     0
## [3,]     1     0     0     0   .     0     0     0     0
## [4,]     0     1     0     0   .     0     0     1     0
## [5,]     0     1     0     0   .     0     0     0     0
## 
## ...
## 
## ,,6
##       [,1]  [,2]  [,3]  [,4] ...  [,8]  [,9] [,10] [,11]
## [1,]     0     2     0     0   .     0     1     0     0
## [2,]     0     1     1     0   .     0     0     1     0
## [3,]     0     1     1     0   .     1     1     0     0
## [4,]     1     0     0     0   .     0     0     0     0
## [5,]     0     0     1     0   .     0     0     0     0

See ?SparseArray_abind for more information and additional examples.

6.4 matrixStats methods

The SparseArray package provides memory-efficient col/row summarization methods for SVT_SparseMatrix objects:

svt7 <- SVT_SparseArray(dim=5:6, dimnames=list(letters[1:5], LETTERS[1:6]))
svt7[c(2, 6, 12:17, 22:30)] <- 101:117

colVars(svt7)
##      A      B      C      D      E      F 
## 2040.2 2080.8 2185.3 3467.0 2443.3    2.5

Note that multidimensional objects are supported:

colVars(svt2)
##      A    B   C
## a   80 1620 NaN
## b  380    0   0
## c  500    0   0
## d 1520 2000 250
colVars(svt2, dims=2)
##        A        B        C 
## 732.6316 857.6316      NaN
colAnyNAs(svt2)
##       A     B     C
## a FALSE FALSE  TRUE
## b FALSE FALSE FALSE
## c FALSE FALSE FALSE
## d FALSE FALSE FALSE
colAnyNAs(svt2, dims=2)
##     A     B     C 
## FALSE FALSE  TRUE

See ?SparseArray_matrixStats for more information and additional examples.

6.5 rowsum() and colsum()

rowsum() and colsum() are supported:

rowsum(svt7, group=c(1:3, 2:1))
##     A   B   C   D   E   F
## 1   0 102 106 107 112 230
## 2 101   0 208 108 220 230
## 3   0   0 104   0 110 115
colsum(svt7, group=c("A", "B", "A", "B", "B", "A"))
##     A   B
## a 113 209
## b 318 217
## c 219 110
## d 221 111
## e 223 112

See ?rowsum_methods for more information and additional examples.

6.6 Matrix multiplication and cross-product

SVT_SparseMatrix objects support matrix multiplication:

svt7 %*% svt4
##   [,1] [,2]
## a  204  220
## b  430  638
## c  434  433
## d  438  437
## e  442  441

as well as crossprod() and tcrossprod():

crossprod(svt4)
##      [,1] [,2]
## [1,]   14    5
## [2,]    5    8

See ?SparseMatrix_mult for more information and additional examples.

7 Other operations

7.1 Generate a random SVT_SparseArray object

Two convenience functions are provided for this:

randomSparseArray(c(5, 6, 2), density=0.5)
## <5 x 6 x 2 SparseArray> of type "double" [nzcount=30 (50%)]:
## ,,1
##        [,1]   [,2]   [,3]   [,4]   [,5]   [,6]
## [1,]  0.000 -0.084 -0.170 -1.300 -2.400  0.000
## [2,]  0.000  0.000 -0.870  0.000  1.700  0.000
## [3,] -1.300  0.000  0.000 -0.690  0.000  1.400
## [4,] -0.300 -0.660 -0.034  0.000  0.000  0.000
## [5,] -1.100  0.750  0.300  0.920  0.130  0.000
## 
## ,,2
##        [,1]   [,2]   [,3]   [,4]   [,5]   [,6]
## [1,]  0.000  1.200  0.000  0.000  0.000  0.000
## [2,]  0.000  0.000  1.200  0.000 -1.600  1.500
## [3,]  0.000  0.000  0.000  0.000  1.300  0.320
## [4,] -0.150  0.000 -1.100  0.000  1.100  0.000
## [5,]  0.000  0.360 -0.110  0.057  0.000 -1.100
poissonSparseArray(c(5, 6, 2), density=0.5)
## <5 x 6 x 2 SparseArray> of type "integer" [nzcount=26 (43%)]:
## ,,1
##      [,1] [,2] [,3] [,4] [,5] [,6]
## [1,]    0    0    0    1    0    2
## [2,]    0    0    0    1    1    0
## [3,]    0    1    0    0    0    0
## [4,]    0    0    0    2    1    1
## [5,]    0    0    0    0    1    0
## 
## ,,2
##      [,1] [,2] [,3] [,4] [,5] [,6]
## [1,]    1    0    1    0    0    1
## [2,]    1    1    0    2    0    0
## [3,]    0    1    2    1    1    2
## [4,]    0    0    0    0    1    0
## [5,]    2    1    1    2    1    0

See ?randomSparseArray for more information and additional examples.

7.2 Read/write a sparse matrix from/to a CSV file

Use writeSparseCSV() to write a sparse matrix to a CSV file:

csv_file <- tempfile()
writeSparseCSV(svt7, csv_file)

Use readSparseCSV() to read the file. This will import the data as an SVT_SparseMatrix object:

readSparseCSV(csv_file)
## <5 x 6 SparseMatrix> of type "integer" [nzcount=17 (57%)]:
##     A   B   C   D   E   F
## a   0 102   0 107   0 113
## b 101   0 103 108 109 114
## c   0   0 104   0 110 115
## d   0   0 105   0 111 116
## e   0   0 106   0 112 117

See ?readSparseCSV for more information and additional examples.

8 Comparison with dgCMatrix objects

8.1 “SVT layout” vs “CSC layout”

The nonzero data of a SVT_SparseArray object is stored in a Sparse Vector Tree. This internal data representation is referred to as the “SVT layout”. It is similar to the “CSC layout” (compressed, sparse, column-oriented format) used by CsparseMatrix derivatives from the Matrix package, like dgCMatrix or lgCMatrix objects, but with the following improvements:

  • The “SVT layout” supports sparse arrays of arbitrary dimensions.

  • With the “SVT layout”, the sparse data can be of any type. Whereas CsparseMatrix derivatives only support sparse data of type "double" or "logical".

  • The “SVT layout” imposes no limit on the number of nonzero array elements that can be stored. With dgCMatrix/lgCMatrix objects, this number must be < 2^31.

  • Overall, the “SVT layout” allows more efficient operations on SVT_SparseArray objects.

See ?SVT_SparseArray for more information about the “SVT layout”.

8.2 Working with a big sparse dataset

The “1.3 Million Brain Cell Dataset” from 10x Genomics is a sparse 2D dataset with more than 2^31 nonzero values. The dataset is stored in an HDF5 file that is available via ExperimentHub (resource id EH1039):

suppressMessages(library(HDF5Array))
suppressMessages(library(ExperimentHub))
hub <- ExperimentHub()
oneM <- TENxMatrix(hub[["EH1039"]], group="mm10")
## see ?TENxBrainData and browseVignettes('TENxBrainData') for documentation
## loading from cache
oneM
## <27998 x 1306127> sparse TENxMatrix object of type "integer":
##                      AAACCTGAGATAGGAG-1 ... TTTGTCATCTGAAAGA-133
## ENSMUSG00000051951                    0   .                    0
## ENSMUSG00000089699                    0   .                    0
## ENSMUSG00000102343                    0   .                    0
## ENSMUSG00000025900                    0   .                    0
## ENSMUSG00000109048                    0   .                    0
##                ...                    .   .                    .
## ENSMUSG00000079808                    0   .                    0
## ENSMUSG00000095041                    1   .                    0
## ENSMUSG00000063897                    0   .                    0
## ENSMUSG00000096730                    0   .                    0
## ENSMUSG00000095742                    0   .                    0

oneM is a TENxMatrix object. This is a particular kind of sparse DelayedArray object where the data is on disk in an HDF5 file. See ?TENxMatrix in the HDF5Array package for more information about TENxMatrix objects.

Note that the object has more than 2^31 nonzero values:

nzcount(oneM)
## [1] 2624828308

The standard way to load the data of a TENxMatrix object (or any DelayedArray derivative) from disk to memory is to simply coerce the object to the desired in-memory representation.

For example, to load the data in an SVT_SparseArray object:

# WARNING: This takes a couple of minutes on a modern laptop, and will
# consume about 25Gb of RAM!
svt <- as(oneM, "SVT_SparseArray")

To load the data in a dgCMatrix object:

# This will fail because 'oneM' has more than 2^31 nonzero values!
as(oneM, "dgCMatrix")

9 Learn more

Please consult the individual man pages in the SparseArray package to learn more about SVT_SparseArray objects and about the package. A good starting point is the man page for SparseArray objects: ?SparseArray

10 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] TENxBrainData_1.25.0        SingleCellExperiment_1.29.0
##  [3] SummarizedExperiment_1.37.0 Biobase_2.67.0             
##  [5] GenomicRanges_1.59.0        GenomeInfoDb_1.43.0        
##  [7] ExperimentHub_2.15.0        AnnotationHub_3.15.0       
##  [9] BiocFileCache_2.15.0        dbplyr_2.5.0               
## [11] HDF5Array_1.35.0            rhdf5_2.51.0               
## [13] DelayedArray_0.33.0         SparseArray_1.7.0          
## [15] S4Arrays_1.7.0              IRanges_2.41.0             
## [17] abind_1.4-8                 S4Vectors_0.45.0           
## [19] MatrixGenerics_1.19.0       matrixStats_1.4.1          
## [21] BiocGenerics_0.53.0         Matrix_1.7-1               
## [23] BiocStyle_2.35.0           
## 
## loaded via a namespace (and not attached):
##  [1] KEGGREST_1.47.0         xfun_0.48               bslib_0.8.0            
##  [4] lattice_0.22-6          rhdf5filters_1.19.0     vctrs_0.6.5            
##  [7] tools_4.5.0             generics_0.1.3          curl_5.2.3             
## [10] tibble_3.2.1            fansi_1.0.6             AnnotationDbi_1.69.0   
## [13] RSQLite_2.3.7           blob_1.2.4              pkgconfig_2.0.3        
## [16] GenomeInfoDbData_1.2.13 lifecycle_1.0.4         compiler_4.5.0         
## [19] Biostrings_2.75.0       htmltools_0.5.8.1       sass_0.4.9             
## [22] yaml_2.3.10             pillar_1.9.0            crayon_1.5.3           
## [25] jquerylib_0.1.4         cachem_1.1.0            mime_0.12              
## [28] tidyselect_1.2.1        digest_0.6.37           purrr_1.0.2            
## [31] dplyr_1.1.4             bookdown_0.41           BiocVersion_3.21.1     
## [34] fastmap_1.2.0           grid_4.5.0              cli_3.6.3              
## [37] magrittr_2.0.3          utf8_1.2.4              withr_3.0.2            
## [40] UCSC.utils_1.3.0        filelock_1.0.3          rappdirs_0.3.3         
## [43] bit64_4.5.2             rmarkdown_2.28          XVector_0.47.0         
## [46] httr_1.4.7              bit_4.5.0               png_0.1-8              
## [49] memoise_2.0.1           evaluate_1.0.1          knitr_1.48             
## [52] rlang_1.1.4             glue_1.8.0              DBI_1.2.3              
## [55] BiocManager_1.30.25     jsonlite_1.8.9          R6_2.5.1               
## [58] Rhdf5lib_1.29.0         zlibbioc_1.53.0