readSparseCSV: Read/write a sparse matrix from/to a CSV file

View source: R/readSparseCSV.R

readSparseCSVR Documentation

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

Description

Read/write a sparse matrix from/to a CSV (comma-separated values) file.

Usage

writeSparseCSV(x, filepath, sep=",", transpose=FALSE, write.zeros=FALSE,
                  chunknrow=250)

readSparseCSV(filepath, sep=",", transpose=FALSE)

Arguments

x

A matrix-like object, typically sparse. IMPORTANT: The object must have rownames and colnames! These will be written to the file.

Another requirement is that the object must be subsettable. More precisely: it must support 2D-style subsetting of the kind x[i, ] and x[ , j] where i and j are integer vectors of valid row and column indices.

filepath

The path (as a single string) to the file where to write the matrix-like object or to read it from. Compressed files are supported.

If "", writeSparseCSV() will write the data to the standard output connection.

Note that filepath can also be a connection.

sep

The field separator character. Values on each line of the file are separated by this character.

transpose

TRUE or FALSE. By default, rows in the matrix-like object correspond to lines in the CSV file. Set transpose to TRUE to transpose the matrix-like object on-the-fly, that is, to have its columns written to or read from the lines in the CSV file.

Note that using transpose=TRUE is semantically equivalent to calling t() on the object before writing it or after reading it, but it will tend to be more efficient. Also it will work even if x does not support t() (not all matrix-like objects are guaranteed to be transposable).

write.zeros

TRUE or FALSE. By default, the zero values in x are not written to the file. Set write.zeros to TRUE to write them.

chunknrow

writeSparseCSV() uses a block-processing strategy to try to speed up things. By default blocks of 250 rows (or columns if transpose=TRUE) are used. In our experience trying to increase this (e.g. to 500 or more) will generally not produce significant benefits while it will increase memory usage, so use carefully.

Value

writeSparseCSV returns an invisible NULL.

readSparseCSV returns a SparseMatrix object of class SVT_SparseMatrix.

See Also

  • SparseArray objects.

  • dgCMatrix objects implemented in the Matrix package.

Examples

## ---------------------------------------------------------------------
## writeSparseCSV()
## ---------------------------------------------------------------------

## Prepare toy matrix 'm0':
rownames0 <- LETTERS[1:6]
colnames0 <- letters[1:4]
m0 <- matrix(0L, nrow=length(rownames0), ncol=length(colnames0),
                 dimnames=list(rownames0, colnames0))
m0[c(1:2, 8, 10, 15:17, 24)] <- (1:8)*10L
m0

## writeSparseCSV():
writeSparseCSV(m0, filepath="", sep="\t")
writeSparseCSV(m0, filepath="", sep="\t", write.zeros=TRUE)
writeSparseCSV(m0, filepath="", sep="\t", transpose=TRUE)

## Note that writeSparseCSV() will automatically (and silently) coerce
## non-integer values to integer by passing them thru as.integer().

## Example where type(x) is "double":
m1 <- m0 * runif(length(m0))
m1
type(m1)
writeSparseCSV(m1, filepath="", sep="\t")

## Example where type(x) is "logical":
writeSparseCSV(m0 != 0, filepath="", sep="\t")

## Example where type(x) is "raw":
m2 <- m0
type(m2) <- "raw"
m2
writeSparseCSV(m2, filepath="", sep="\t")

## ---------------------------------------------------------------------
## readSparseCSV()
## ---------------------------------------------------------------------

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

svt1 <- readSparseCSV(csv_file)
svt1

svt2 <- readSparseCSV(csv_file, transpose=TRUE)
svt2

## If you need the sparse data as a dgCMatrix object, just coerce the
## returned object:
as(svt1, "dgCMatrix")
as(svt2, "dgCMatrix")

## Sanity checks:
stopifnot(identical(m0, as.matrix(svt1)))
stopifnot(identical(t(m0), as.matrix(svt2)))

Bioconductor/SparseArray documentation built on Nov. 14, 2024, 6:22 p.m.