matmult | R Documentation |
Multithreaded <matrix, matrix> multiplications ('%*%', 'crossprod', and 'tcrossprod') and <matrix, vector> multiplications ('%*%'), for <sparse, dense> matrix combinations and <sparse, vector> combinations (See signatures for supported combinations).
Objects from the 'float' package are also supported for some combinations.
## S4 method for signature 'matrix,CsparseMatrix'
x %*% y
## S4 method for signature 'float32,CsparseMatrix'
x %*% y
## S4 method for signature 'matrix,RsparseMatrix'
tcrossprod(x, y)
## S4 method for signature 'float32,RsparseMatrix'
tcrossprod(x, y)
## S4 method for signature 'matrix,CsparseMatrix'
crossprod(x, y)
## S4 method for signature 'float32,CsparseMatrix'
crossprod(x, y)
## S4 method for signature 'RsparseMatrix,matrix'
tcrossprod(x, y)
## S4 method for signature 'RsparseMatrix,matrix'
x %*% y
## S4 method for signature 'RsparseMatrix,float32'
x %*% y
## S4 method for signature 'RsparseMatrix,float32'
tcrossprod(x, y)
## S4 method for signature 'RsparseMatrix,numeric'
x %*% y
## S4 method for signature 'RsparseMatrix,logical'
x %*% y
## S4 method for signature 'RsparseMatrix,integer'
x %*% y
## S4 method for signature 'RsparseMatrix,sparseVector'
x %*% y
x , y |
dense ( |
Will try to use the maximum available number of threads for the computations when appropriate. The number of threads can be controlled through the package options (e.g. 'options("MatrixExtra.nthreads" = 1)' - see MatrixExtra-options) and will be set to 1 after running restore_old_matrix_behavior.
Be aware that sparse-dense matrix multiplications might suffer from reduced numerical precision, especially when using objects of type 'float32' (from the 'float' package).
Internally, these functions use BLAS level-1 routines, so their speed might depend on the BLAS backend being used (e.g. MKL, OpenBLAS) - that means: they might be quite slow on a default install of R for Windows (see this link for a tutorial about getting OpenBLAS in R for Windows).
Doing computations in float32 precision depends on the package float, and as such comes with some caveats:
On Windows, if installing 'float' from CRAN, it will use very unoptimized routines which will likely result in a slowdown compared to using regular double (numeric) type. Getting it to use an optimized BLAS library is not as simple as substituting the Rblas DLL - see the package's README for details.
On macOS, it will use static linking for 'float', thus if changing the BLAS library used by R, it will not change the float32 functions, and getting good performance out of it might require compiling it from source with '-march=native' flag.
When multiplying a sparse matrix by a sparse vector, their indices will be sorted in-place (see sort_sparse_indices).
In order to match exactly with base R's behaviors, when passing vectors to these operators, will assume their shape as follows:
MatMult(Matrix, vector): column vector if the matrix has more than one column or is empty, row vector if the matrix has only one column.
MatMult(vector, Matrix): row vector if the matrix has more than one row, column vector if the matrix has only one row.
MatMul(vector, vector): LHS is a row vector, RHS is a column vector.
crossprod(Matrix, vector): column vector if the matrix has more than one row, row vector if the matrix has only one row.
crossprod(vector, Matrix): column vector.
crossprod(vector, vector): column vector.
tcrossprod(Matrix, vector): row vector if the matrix has only one row, column vector if the matrix has only one column, and will throw an error otherwise.
tcrossprod(vector, Matrix): row vector if the matrix has more than one column, column vector if the matrix has only one column.
tcrossprod(vector, vector): column vector.
In general, the output returned by these functions will be a dense matrix from base R, or a dense matrix from 'float' when one of the inputs is also from the 'float' package, with the following exceptions:
MatMult(RsparseMatrix[n,1], vector) -> 'dgRMatrix'.
MatMult(RsparseMatrix[n,1], sparseVector) -> 'dgCMatrix'.
MatMult(float32[n], CsparseMatrix[1,m]) -> 'dgCMatrix'.
tcrossprod(float32[n], RsparseMatrix[m,1]) -> 'dgCMatrix'.
A dense matrix
object in most cases, with some exceptions which might
come in sparse format (see the 'Details' section).
library(Matrix)
library(MatrixExtra)
### To use all available threads (default)
options("MatrixExtra.nthreads" = parallel::detectCores())
### Example will run with only 1 thread (CRAN policy)
options("MatrixExtra.nthreads" = 1)
## Generate random matrices
set.seed(1)
A <- rsparsematrix(5,4,.5)
B <- rsparsematrix(4,3,.5)
## Now multiply in some supported combinations
as.matrix(A) %*% as.csc.matrix(B)
as.csr.matrix(A) %*% as.matrix(B)
crossprod(as.matrix(B), as.csc.matrix(B))
tcrossprod(as.csr.matrix(A), as.matrix(A))
### Restore the number of threads
options("MatrixExtra.nthreads" = parallel::detectCores())
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.