Description of functions to perform matrix operations, algebra and some basic statistical models using Matrix (CRAN) and DelayedMatrix (Bioconductor) objects.
BigDataStatMeth 0.99.14
This package implements several matrix operations using Matrix
and DelayMatrix
objects as well as HDF5 data files. Some basic algebra operations that can also be computed that are useful to implement statistical analyses using standard methodologies such as principal component analyses (PCA) or least squares estimation. The package also contains specific statistical methods mainly used in omic
data analysis such as lasso regression. All procedures referred to HDF5 can be found in BigDataStatMeth_hdf5 vignette.
The package requires other packages from CRAN
and Bioconductor
to be installed.
CRAN
: Matrix
, RcppEigen
and RSpectra
.Bioconductor
: beachmat
, DelayedArray
As the package can also deal with hdf5 files [See Vignette BigDataStatMeth_hdf5], these other packages are required: HDF5Array
, rhdf5
. The user can execute this code to install the required packages:
# Install BiocManager (if not previously installed)
install.packages("BiocManager")
# Install required packages
BiocManager::install(c("Matrix", "RcppEigen", "RSpectra",
"beachmat", "DelayedArray",
"HDF5Array", "rhdf5"))
Our package needs to be installed from source code. In such cases, a collection of software (e.g. C, C++, Fortran, etc.) are required, mainly for Windows users. These programs can be installed using Rtools.
Once required packages and Rtools are installed, BigDataStatMeth
package can be installed from our GitHub repository as follows:
# Install devtools and load library (if not previously installed)
install.packages("devtools")
library(devtools)
# Install BigDataStatMeth
install_github("isglobal-brge/BigDataStatMeth")
First, let us start by loading the required packages to describe the main capabilities of the package
library(Matrix)
library(DelayedArray)
Warning: package 'MatrixGenerics' was built under R version 4.1.1
Warning: package 'S4Vectors' was built under R version 4.1.1
library(BigDataStatMeth)
library(ggplot2)
This packages is also required to reproduce this vignette
library(microbenchmark)
DelayedMatrix
A DelayedMatrix
is a type of DelayedArray
which is like an ordinary array in R, but allows for the data to be in-memory, on-disk in a file, or even hosted on a remote server.
The DelayedArray framework supports some on-disk backends for genomic data. We can deal with HDF5 files via rhdf5 and the HDF5Array package, GenomicDataStorage (GDS) via the gdsfmt and GDSArray package, and Variant Calling Format (VCF) files via the VCFArray package. However, this can be extended to support other on-disk backends. In theory, it should be possible to implement a DelayedArray backend for any file format that has the capability to store array data with fast random access.
Having data in such format allows to perform very fast operations for a given matrix. This has been implemented in the DelayedMatrixStats package. Several examples of how to deal with DelayedArray
can be found in this workshop.
In this section, different products of matrices and vectors are introduced. The methods implement different strategies including block multiplication algorithms and the use of parallel implementations.
A block matrix or a partitioned matrix is a matrix that is interpreted as having been broken into sections called blocks or submatrices. Intuitively, a block matrix can be visualized as the original matrix with a collection of horizontal and vertical lines, which break it up, or partition it, into a collection of smaller matrices. the implementation has been made from the adaptation of the Fox algorithm [1].
\[A*B=\begin{pmatrix} {A}_{11}&{A}_{12} \\ {A}_{21}&{A}_{22}\end{pmatrix}*\begin{pmatrix}{B}_{11}&B_{12}\\B_{21}&B_{22}\end{pmatrix}=\begin{pmatrix}{C}_{11}&C_{12}\\C_{21}&C_{22}\end{pmatrix}=C\]
Let us simulate to set of matrices to illustrate the use of the function accross the entire documment. First, we simulate a simple case with to matrices A and B with dimensions 500x500 and 500x750, respectively. Second, another example with dimensions 1000x10000 are use to evaluate the performance in large matrices. The examples with big datasets will be illustrated using real data belonging to different omic settings. We can simulate to matrices with the desired dimensions by
# Define small matrix A and B
set.seed(123456)
n <- 500
p <- 750
A <- matrix(rnorm(n*n), nrow=n, ncol=n)
B <- matrix(rnorm(n*p), nrow=n, ncol=p)
# Define big matrix Abig and Bbig
n <- 1000
p <- 10000
Abig <- matrix(rnorm(n*n), nrow=n, ncol=n)
Bbig <- matrix(rnorm(n*p), nrow=n, ncol=p)
They can be converted into DelayedMatrix
object by simply
DA <- DelayedArray(A)
DB <- DelayedArray(B)
Matrix multiplication using block matrices is implemented in the bdblockmult()
function. It only requires to setting the argument block_size
, by default block_size = 128
. An optimal block size is important to better performance but it is difficult to asses what is the optimum block size for each matrix.
# Use 10x10 blocks
AxB <- bdblockmult(A, B, block_size = 10)
# Use 256x256 blocks
AxBDelay <- bdblockmult(DA, DB, block_size = 256 )
As expected the results obtained using this procedure are the correct ones
all.equal(AxBDelay,A%*%B)
[1] TRUE
all.equal(AxB, AxBDelay)
[1] TRUE
Note that when the argument block_size
is larger than any of the dimensions of matrix A or B the blocks_size
is set to min(cols(A), rows(A), cols(B), rows(B))
.
The process can be speed up by making computations in parallel using paral=TRUE
an optional parameter threads can be used to indicate the number of threads
to launch simultaneously, if threads=NULL
the function takes the available threads - 1, leaving one available for user.
AxB <- bdblockmult(A, B, block_size = 10, paral = TRUE)
AxBDelay <- bdblockmult(DA, DB, block_size = 10, paral = TRUE )
AxBDelayT <- bdblockmult(DA, DB, block_size = 10, paral = TRUE, threads = 3 )
all.equal(AxBDelay,A%*%B)
[1] TRUE
all.equal(AxB, AxBDelay)
[1] TRUE
all.equal(AxBDelay, AxBDelayT)
[1] TRUE
To work with big matrices bdblockmult()
can save matrices in hdf5 file format in order to be able to operate with them in blocks and not overload the memory, by default are considered large matrices if the number of rows or columns is greater than 5000, but it can be changed with bigmatrix
argument. Or we can also force to execute the matrix multiplication with data on memory setting the parameter onmemory = TRUE
DAbig <- DelayedArray(Abig)
DBbig <- DelayedArray(Abig)
# We want to force it to run in memory
AxBNOBig <- bdblockmult(Abig, Bbig, onmemory = TRUE)
# Run matrix multiplication with data on memory using submatrices of 256x256
AxBBig3000 <- bdblockmult(DAbig, DBbig, block_size = 256 , onmemory = TRUE)
Here, we compare the performance of the block method with the different options.
bench1 <- microbenchmark(
# Parallel block size = 128
Paral128Mem = bdblockmult(Abig, Bbig, paral = TRUE),
# On disk parallel block size = 256
Paral256Disk = bdblockmult(Abig, Bbig, block_size=256, paral=TRUE),
Paral256Mem = bdblockmult(Abig, Bbig, block_size=256,
paral=TRUE, onmemory=TRUE),
Paral1024Mem = bdblockmult(Abig, Bbig, block_size=1024,
paral=TRUE, onmemory=TRUE), times = 3 )
bench1
Unit: seconds
expr min lq mean median uq max neval cld
Paral128Mem 8.134192 8.367881 8.772597 8.601570 9.091800 9.582029 3 c
Paral256Disk 4.983364 5.034852 5.220670 5.086340 5.339322 5.592305 3 b
Paral256Mem 1.755980 1.757631 1.786949 1.759281 1.802434 1.845586 3 a
Paral1024Mem 1.574793 1.589548 1.601708 1.604303 1.615166 1.626029 3 a
We can observe that the shortest execution time is achieved with block_size = 256
.
The same information is depicted in the next Figure which illustrates the comparison between the different assessed methods
ggplot2::autoplot(bench1)
A sparse matrix or sparse array is a matrix in which most of the elements are zero. BigDataStatMeth
allows to perform matrix multiplication with sparse matrices using the function bdblockmult_sparse
. It is necessary that at least one of the the two matrices is defined as sparse in R using dgCMatrix
class. An example of a sparse matrix could be :
\[\begin{equation} \begin{bmatrix} a_{11}&a_{12} & 0 & \cdots & \cdots & \cdots & \cdots & 0 \\ a_{21} & a_{22} & a_{23} & \ddots & && & \vdots \\ 0 & a_{32} & a_{33} & a_{34} & \ddots & & & \vdots \\ \vdots & \ddots & \ddots & \ddots & \ddots & \ddots & & \vdots \\ \vdots & & \ddots & \ddots & \ddots & \ddots & \ddots& \vdots\\ \vdots & & & & \ddots & a_{n-1,n-2} & a_{n-1,n-1} & a_{n-1,n}\\ 0 & \cdots & \cdots & \cdots & \cdots & 0 & a_{n,n-1} & a_{n,n} \\ \end{bmatrix} \end{equation}\]
k <- 1e3
# Generate 2 sparse matrix x_sparse and y_sparse
set.seed(1)
x_sparse <- sparseMatrix(
i = sample(x = k, size = k),
j = sample(x = k, size = k),
x = rnorm(n = k)
)
set.seed(2)
y_sparse <- sparseMatrix(
i = sample(x = k, size = k),
j = sample(x = k, size = k),
x = rnorm(n = k)
)
d <- bdblockmult_sparse(x_sparse, y_sparse)
f <- x_sparse%*%y_sparse
all.equal(d,f)
[1] TRUE
Here, we compare the performace using sparse matrix multiplication with sparse matrices and the blockmult multiplication with the same matrices not declared as sparce.
res <- microbenchmark(
sparse_mult = bdblockmult_sparse(x_sparse, y_sparse),
matrix_mult = bdblockmult(as.matrix(x_sparse), as.matrix(y_sparse)),
RSparse_mult = x_sparse%*% y_sparse,
times = 3 )
res
Unit: microseconds
expr min lq mean median uq
sparse_mult 53.144 56.031 69.1790 58.918 77.1965
matrix_mult 194844.093 195859.694 196324.1323 196875.296 197064.1520
RSparse_mult 115.666 120.807 155.1513 125.948 174.8940
max neval cld
95.475 3 a
197253.008 3 b
223.840 3 a
The same information is depicted in the next Figure which illustrates the comparison between the different assessed methods
ggplot2::autoplot(res)
We can observe the huge difference in execution time
To perform a cross-product and tcrossproduct with BigDataStatMeth we use the same function bdcrossprod()
setting parameter transposed
to TRUE or FALSE. Like other functions implemented in BigDataStatMeth, we can work with DelayedArray or R Objects. Setting parameter transposed for Crossproduct and tCrossProduct
if transposed = FALSE
(default), we perform a Crossproduct \[\left( C \right) ={ \left( A \right) }^{ t }\left( A \right) \]
if transposed = TRUE
: we perform a transposed-Crossproduct \[\left( C \right) = \left( A \right){ \left( A \right) }^{ t } \]
Here we shown some examples using bdcrossprod
function
n <- 500
m <- 250
A <- matrix(rnorm(n*m), nrow=n, ncol=m)
DA <- DelayedArray(A)
# Cross Product of a standard R matrix
cpA <- bdCrossprod(A)
# Result with DelayedArray data type
cpDA <- bdCrossprod(DA)
We obtain the expected values computed using crossprod
R function
all.equal(cpDA, crossprod(A))
[1] TRUE
you may also set transposed=TRUE
# Transposed Cross Product R matrices
tcpA <- bdtCrossprod(A)
# With DelayeArray data types
tcpDA <- bdtCrossprod(DA)
We obtain the expected values computed using tcrossprod
function
all.equal(tcpDA, tcrossprod(A))
[1] TRUE
We can show that the implemented version really improves the R implementation computational speed.
res <- microbenchmark(
bdcrossp_tr = bdtCrossprod(A),
rcrossp_tr = tcrossprod(A),
bdcrossp = bdCrossprod(A),
rcrossp = crossprod(A),
times = 3)
res
Unit: milliseconds
expr min lq mean median uq max neval
bdcrossp_tr 6.921249 7.221568 7.512567 7.521887 7.808226 8.094566 3
rcrossp_tr 21.528390 21.543555 21.577063 21.558720 21.601400 21.644080 3
bdcrossp 2.942669 2.968552 3.152647 2.994436 3.257636 3.520836 3
rcrossp 18.346368 18.348487 18.392689 18.350607 18.415850 18.481092 3
cld
b
d
a
c
ggplot2::autoplot(res)
You can perform a weighted cross-product \(C = X^ t w X\) with bdwcrossprod()
given a matrix X as argument and a vector or matrix of weights, w.
n <- 250
X <- matrix(rnorm(n*n), nrow=n, ncol=n)
u <- runif(n)
w <- u * (1 - u)
DX <- DelayedArray(X)
Dw <- DelayedArray(as.matrix(w))
wcpX <- bdwproduct(X, w,"xwxt")
wcpDX <- bdwproduct(DX, Dw,"xwxt") # with DelayedArray
wcpDX[1:5,1:5]
[,1] [,2] [,3] [,4] [,5]
[1,] 44.1890362 1.032160 0.7061483 -1.8762456 -1.2572584
[2,] 1.0321595 43.919531 -2.2856272 -1.3729013 1.1013293
[3,] 0.7061483 -2.285627 46.4180971 -3.6621389 2.1444686
[4,] -1.8762456 -1.372901 -3.6621389 39.1658810 -0.5414553
[5,] -1.2572584 1.101329 2.1444686 -0.5414553 39.1266511
Those are the expected values as it is indicated by executing:
all.equal( wcpDX, X%*%diag(w)%*%t(X) )
[1] TRUE
With argument transposed=TRUE
, we can perform a transposed weighted cross-product \(C = A w A^t\)
wtcpX <- bdwproduct(X, w,"xtwx")
wtcpDX <- bdwproduct(DX, Dw,"xtwx") # with DelayedArray
wtcpDX[1:5,1:5]
[,1] [,2] [,3] [,4] [,5]
[1,] 42.10759651 -0.05701354 -0.4767347 -6.3149197 -0.4659907
[2,] -0.05701354 50.97505787 -3.9501323 -0.7959214 -0.0952895
[3,] -0.47673475 -3.95013234 43.4955384 -0.5409623 4.2215983
[4,] -6.31491965 -0.79592141 -0.5409623 43.0341187 -2.0391817
[5,] -0.46599066 -0.09528950 4.2215983 -2.0391817 43.7745996
Those are the expected values as it is indicated by executing:
all.equal(wtcpDX, t(X)%*%diag(w)%*%X)
[1] TRUE
The Cholesky factorization is widely used for solving a system of linear equations whose coefficient matrix is symmetric and positive definite.
\[A = LL^t = U^tU\]
where \(L\) is a lower triangular matrix and U is an upper triangular matrix. To get the Inverse Cholesky we can use the function bdInvCholesky()
. Let us start by illustrating how to do this computations in a simulated data:
# Generate a positive definite matrix
Posdef <- function (n, ev = runif(n, 0, 10))
{
Z <- matrix(ncol=n, rnorm(n^2))
decomp <- qr(Z)
Q <- qr.Q(decomp)
R <- qr.R(decomp)
d <- diag(R)
ph <- d / abs(d)
O <- Q %*% diag(ph)
Z <- t(O) %*% diag(ev) %*% O
return(Z)
}
A <- Posdef(n = 500, ev = 1:500)
DA <- DelayedArray(A)
invchol <- bdInvCholesky(A)
Dinvchol <- bdInvCholesky(DA)
round(invchol[1:5,1:5],8)
[,1] [,2] [,3] [,4] [,5]
[1,] 0.01055505 0.00087913 -0.00090136 0.00442774 0.00074769
[2,] 0.00087913 0.01467257 -0.00442116 0.00714249 0.00001366
[3,] -0.00090136 -0.00442116 0.01744738 -0.00044069 -0.00352095
[4,] 0.00442774 0.00714249 -0.00044069 0.03232484 -0.00420970
[5,] 0.00074769 0.00001366 -0.00352095 -0.00420970 0.01388255
We can check whether this function returns the expected values obtained with the standard R function solve
:
all.equal(Dinvchol, solve(A))
[1] TRUE
We can show that the implemented version really improves the R implementation computational speed.
res <- microbenchmark(invchol = bdInvCholesky(A),
invcholDelayedA = bdInvCholesky(DA),
invcholR = solve(A),
times = 3)
res
Unit: milliseconds
expr min lq mean median uq max neval
invchol 36.75295 37.30357 39.36938 37.85419 40.67760 43.50101 3
invcholDelayedA 39.36129 39.91267 41.05822 40.46406 41.90669 43.34932 3
invcholR 86.86038 87.38683 88.01469 87.91329 88.59184 89.27040 3
cld
a
a
b
ggplot2::autoplot(res)
The Moore-Penrose pseudoinverse is a direct application of the SVD. The inverse of a matrix A can be used to solve the equation \({Ax}={b}\). But in the case where the set of equations have 0 or many solutions the inverse cannot be found and the equation cannot be solved. The following formula is used to find the pseudoinverse:
\[{A}^+= {VD}^+{U}^T\]
In BigDataStatMeth we implemented Moore-Penrose pseudoinverse in function bdpseudoinv
, now we get the pseudoinverse from a simulated data
m <- 1300
n <- 1200
A <- matrix(rnorm(n*m), nrow=n, ncol=m)
pseudoinv <- bdpseudoinv(A)
zapsmall(pseudoinv)[1:5,1:5]
[,1] [,2] [,3] [,4] [,5]
[1,] -0.002881082 0.003254617 0.006027423 0.001720531 -0.000531196
[2,] 0.001387901 0.002703614 -0.003206459 -0.000507613 -0.000079770
[3,] 0.002056213 -0.000729940 -0.001762621 0.004061003 0.001203367
[4,] -0.004359027 0.001531990 -0.000333149 -0.000417425 -0.001130834
[5,] 0.000029255 -0.000038933 0.001082334 -0.000678160 0.002372722
QR decomposition, also known as a QR factorization or QU factorization is a decomposition of a matrix A into a product : \[A = QR\] of an orthogonal matrix Q and an upper triangular matrix R. QR decomposition is often used to solve the linear least squares problems.
In BigDataStatMeth we implemented QR decomposition in function bdQR
, the QR decomposition can be performed in R objects or in DelayedArray objects. To show how to use this function we performa a QR decomposition from the previous simulated data in matrix A.
# Get Delayed array from A
DA <- DelayedArray(A)
QR_A <- bdQR(A)
QR_DA <- bdQR(DA)
QR_R <- qr(A)
# Show results for Q
zapsmall(QR_A$Q[1:5,1:5])
[,1] [,2] [,3] [,4] [,5]
[1,] -0.04130003 -0.01870761 -0.02889834 0.00313758 -0.03594914
[2,] -0.01186500 -0.00624750 -0.02207752 0.01768278 0.02601456
[3,] 0.01540677 0.01364885 -0.02911809 0.00598247 -0.04110039
[4,] -0.03217148 0.01177411 -0.00549290 -0.04301058 -0.02591033
[5,] -0.02398464 -0.00726357 -0.00909492 -0.02443984 -0.04444592
# Show results for R
zapsmall(QR_A$R[1:5,1:5])
[,1] [,2] [,3] [,4] [,5]
[1,] 34.967 -0.61199 -1.42200 -2.19133 -0.44725
[2,] 0.000 -34.95367 -0.39611 -0.13508 1.19882
[3,] 0.000 0.00000 -34.09566 0.18723 -0.29108
[4,] 0.000 0.00000 0.00000 -34.24813 -1.16650
[5,] 0.000 0.00000 0.00000 0.00000 -34.10602
# Test Q
all.equal(QR_A$Q, qr.Q(QR_R), check.attributes=FALSE)
[1] TRUE
In BigDataStatMeth we implemented the function bdSolve
that computes the solution to a real system of linear equations
\[A*X = B\] where A is an N-by-N matrix and X and B are N-by-K matrices.
Here we solve a matrix equation with a squared matrix A (1000 by 1000) and B matrix (1000 by 2)
# Simulate data
m <- 1000
n <- 1000
A <- matrix(runif(n*m), nrow = n, ncol = m)
B <- matrix(runif(n*2), nrow = n)
# Solve matrix equation
X <- bdSolve(A, B)
# Show results
X[1:5,]
[,1] [,2]
[1,] -0.7065915 -0.4135584
[2,] -1.0728821 -0.4345175
[3,] 2.1275535 0.1552309
[4,] 1.3388455 0.1721595
[5,] -0.5465846 0.1982003
Now we check results multiplying matrix A by results X, if all is correct \(A*X = B\)
testB <- bdblockmult(A,X)
B[1:5,]
[,1] [,2]
[1,] 0.1126148 0.8005922
[2,] 0.3915829 0.1997295
[3,] 0.7461557 0.1041076
[4,] 0.4530352 0.5773070
[5,] 0.3525574 0.8557561
testB[1:5,]
[,1] [,2]
[1,] 0.1126148 0.8005922
[2,] 0.3915829 0.1997295
[3,] 0.7461557 0.1041076
[4,] 0.4530352 0.5773070
[5,] 0.3525574 0.8557561
all.equal(B, testB)
[1] TRUE
The SVD of an \(m \times n\) real or complex matrix \(A\) is a factorization of the form:
\[U\Sigma { V }^{ T }\]
where :
-\(U\) is a \(m \times m\) real or complex unitary matrix -\(\Sigma\) is a \(m \times n\) rectangular diagonal matrix with non-negative real numbers on the diagonal -\(V\) is a \(n \times n\) real or complex unitary matrix.
Notice that:
He have implemented the SVD for R matrices and Delayed Array objects in the function bdSVD()
. The method, so far, only allows SVD of real matrices \(A\). This code illustrate how to perform such computations:
# Matrix simulation
set.seed(413)
n <- 500
A <- matrix(rnorm(n*n), nrow=n, ncol=n)
# Get a Delayed Array object
DA <- DelayedArray(A)
# SVD
bsvd <- bdSVD(A)
Dbsvd <- bdSVD(DA)
# Singular values, and right and left singular vectors
bsvd$d[1:5]
[1] 44.27037 43.95609 43.31037 43.01980 42.81728
bsvd$u[1:5,1:5]
[,1] [,2] [,3] [,4] [,5]
[1,] 0.028712445 0.10169472 0.019783036 0.032804622 -0.048915416
[2,] 0.001193891 0.02296487 0.009024529 0.059683600 0.027888834
[3,] -0.008878497 -0.02414059 -0.006035909 -0.004867886 0.024570252
[4,] 0.037323916 -0.04675198 0.035928628 -0.021097593 -0.073485962
[5,] -0.054073083 0.04382626 -0.009027050 -0.002269443 -0.003071978
bsvd$v[1:5,1:5]
[,1] [,2] [,3] [,4] [,5]
[1,] -0.07105335 0.04554495 -0.03850442 0.07270812 0.078124410
[2,] 0.05989691 -0.02957792 -0.07508850 -0.01171385 0.001581823
[3,] -0.00200688 -0.04179449 0.02798939 -0.05686766 -0.073573833
[4,] -0.01284403 -0.05024617 -0.09038993 -0.05725360 0.023632656
[5,] -0.01751359 0.07320492 -0.02120536 -0.02397367 -0.067643671
We get the expected results obtained with standard R functions:
all.equal( sqrt( svd( tcrossprod( scale(A) ) )$d[1:10] ), bsvd$d[1:10] )
[1] TRUE
all.equal( sqrt( svd( tcrossprod( scale(A) ) )$d[1:10] ), Dbsvd$d[1:10] )
[1] TRUE
you get the \(\sigma_i\), \(U\) and \(V\) of normalized matrix \(A\), if you want to perform the SVD from not normalized matrix \(A\) then you have to set the parameter bcenter = false
and bscale = false
.
bsvd <- bdSVD(A, bcenter = FALSE, bscale = FALSE)
Dbsvd <- bdSVD(DA, bcenter = FALSE, bscale = FALSE)
bsvd$d[1:5]
[1] 44.44007 43.89640 43.38384 43.23563 42.82658
bsvd$u[1:5,1:5]
[,1] [,2] [,3] [,4] [,5]
[1,] 3.075740e-02 0.09861737 -0.026705796 0.001569226 -0.04384078
[2,] 6.480471e-05 0.04151230 -0.049215808 0.054763599 -0.01436232
[3,] -1.857597e-02 -0.02463540 -0.001246071 0.015838141 0.05725309
[4,] 3.555600e-02 -0.05715730 -0.009596880 -0.040247700 -0.03888504
[5,] -5.290501e-02 0.04784627 -0.006624289 0.010168748 0.01934154
bsvd$v[1:5,1:5]
[,1] [,2] [,3] [,4] [,5]
[1,] -0.0677291931 0.03646765 -0.072013679 -0.06460296 0.03834539
[2,] 0.0558086446 -0.02816709 -0.067619533 0.02478961 0.01316491
[3,] 0.0024305441 -0.04679000 0.064115513 0.04134591 -0.03519768
[4,] -0.0007861777 -0.06602592 -0.062833191 0.06307396 0.03651487
[5,] -0.0165396063 0.07316212 0.002600021 0.07030477 -0.07916813
all.equal( sqrt(svd(tcrossprod(A))$d[1:10]), bsvd$d[1:10] )
[1] TRUE
all.equal( sqrt(svd(tcrossprod(A))$d[1:10]), Dbsvd$d[1:10] )
[1] TRUE
A method developed by M. A. Iwen and B. W. Ong uses a distributed and incremental SVD algorithm that is useful for agglomerative data analysis on large networks. The algorithm calculates the singular values and left singular vectors of a matrix A, by first, partitioning it by columns. This creates a set of submatrices of A with the same number of rows, but only some of its columns. After that, the SVD of each of the submatrices is computed. The final step consists of combining the results obtained by merging them again and computing the SVD of the resulting matrix.
This approach is only implemented using HDF5 files. This method is implemented in bdSVD_hdf5
function, this function works directly on hdf5 data format, loading in memory only the data to perform calculations and saving the results again in the hdf5 file for later use. The user is referred to read section 7.3.1 from this vignette.
DelayedMatrix
objects to implement basic statistical methodsLet us illustrate how to perform a PCA using miRNA data obtained from TCGA corresponding to 3 different tumors: melanoma (ME), leukemia (LEU) and centran nervous system (CNS). Data is available at BigDataStatMeth
and can be loaded by simply:
data(miRNA)
data(cancer)
dim(miRNA)
[1] 21 537
We observe that there are a total of 21 individuals and 537 miRNAs. The vector cancer
contains the type of tumor of each individual. For each type we have:
table(cancer)
cancer
CNS LE ME
6 6 9
Now, the typical principal component analysis on the samples
can be run on the miRNA
matrix since it has miRNAs in columns and individuals in rows
pc <- prcomp(miRNA)
We can plot the two first components with:
plot(pc$x[, 1], pc$x[, 2],
main = "miRNA data on tumor samples",
xlab = "PC1", ylab = "PC2", type="n")
abline(h=0, v=0, lty=2)
points(pc$x[, 1], pc$x[, 2], col = cancer,
pch=16, cex=1.2)
legend("topright", levels(cancer), pch=16, col=1:3)
The same analysis can be performed using SVD decomposition and having miRNAs as a DelayedMatrix
object. The PCA is equivalent to performing the SVD on the centered data, where the centering occurs on the columns. In that case the function bdSVD
requires to set the argument bcenter
and bscale
equal to TRUE, the dafault values. Notice that sweep()
and colMeans()
functions can be applied to a DelayedMAtrix
object since this method is implemented for that type of objects in the DelayedArray
package.
miRNAD <- DelayedArray(miRNA)
miRNAD.c <- DelayedArray::sweep(miRNAD, 2,
DelayedArray::colMeans(miRNAD), "-")
svd.da <- bdSVD(miRNAD.c, bcenter = FALSE, bscale = FALSE)
The PCA plot for the two principal components can then be be obtained by:
plot(svd.da$u[, 1], svd.da$u[, 2],
main = "miRNA data on tumor samples",
xlab = "PC1", ylab = "PC2", type="n")
abline(h=0, v=0, lty=2)
points(svd.da$u[, 1], svd.da$u[, 2], col = cancer,
pch=16, cex=1.2)
legend("topright", levels(cancer), pch=16, col=1:3)
We can observe that both figures are equal irrespective to a sign change of second component (that can happen in SVD).
sessionInfo()
R version 4.1.0 Patched (2021-06-26 r80567)
Platform: x86_64-apple-darwin17.0 (64-bit)
Running under: macOS Mojave 10.14.6
Matrix products: default
BLAS: /Library/Frameworks/R.framework/Versions/4.1/Resources/lib/libRblas.0.dylib
LAPACK: /Library/Frameworks/R.framework/Versions/4.1/Resources/lib/libRlapack.dylib
locale:
[1] C/ca_ES.UTF-8/ca_ES.UTF-8/C/ca_ES.UTF-8/ca_ES.UTF-8
attached base packages:
[1] parallel stats4 stats graphics grDevices utils datasets
[8] methods base
other attached packages:
[1] microbenchmark_1.4-7 ggplot2_3.3.5 DelayedArray_0.18.0
[4] IRanges_2.26.0 S4Vectors_0.30.2 MatrixGenerics_1.4.3
[7] matrixStats_0.61.0 BiocGenerics_0.38.0 Matrix_1.3-4
[10] rmarkdown_2.11 BigDataStatMeth_0.99.14 rhdf5_2.36.0
[13] knitr_1.36 BiocStyle_2.20.2
loaded via a namespace (and not attached):
[1] Rcpp_1.0.7 mvtnorm_1.1-3 lattice_0.20-45
[4] zoo_1.8-9 assertthat_0.2.1 digest_0.6.28
[7] utf8_1.2.2 R6_2.5.1 evaluate_0.14
[10] highr_0.9 pillar_1.6.4 rlang_0.4.12
[13] multcomp_1.4-17 data.table_1.14.2 jquerylib_0.1.4
[16] magick_2.7.3 splines_4.1.0 stringr_1.4.0
[19] RCurl_1.98-1.5 munsell_0.5.0 beachmat_2.8.1
[22] compiler_4.1.0 xfun_0.28.7 pkgconfig_2.0.3
[25] htmltools_0.5.2 tidyselect_1.1.1 tibble_3.1.6
[28] bookdown_0.24 codetools_0.2-18 fansi_0.5.0
[31] crayon_1.4.2 dplyr_1.0.7 withr_2.4.2
[34] MASS_7.3-54 bitops_1.0-7 rhdf5filters_1.4.0
[37] grid_4.1.0 jsonlite_1.7.2 gtable_0.3.0
[40] lifecycle_1.0.1 DBI_1.1.1 magrittr_2.0.1
[43] scales_1.1.1 RcppParallel_5.1.4 stringi_1.7.5
[46] farver_2.1.0 bslib_0.3.1 ellipsis_0.3.2
[49] vctrs_0.3.8 generics_0.1.1 sandwich_3.0-1
[52] TH.data_1.1-0 Rhdf5lib_1.14.2 tools_4.1.0
[55] glue_1.5.0 purrr_0.3.4 survival_3.2-13
[58] fastmap_1.1.0 yaml_2.2.1 colorspace_2.0-2
[61] BiocManager_1.30.16 sass_0.4.0