write10xCounts {DropletUtils}R Documentation

Write count data in the 10x format

Description

Create a directory containing the count matrix and cell/gene annotation from a sparse matrix of UMI counts, in the format produced by the CellRanger software suite.

Usage

write10xCounts(path, x, barcodes=colnames(x), gene.id=rownames(x), 
    gene.symbol=gene.id, overwrite=FALSE)

Arguments

x

A sparse numeric matrix of UMI counts.

path

A string containing the path to the output directory.

barcodes

A character vector of cell barcodes, one per column of x.

gene.id

A character vector of gene identifiers, one per row of x.

gene.symbol

A character vector of gene symbols, one per row of x.

overwrite

A logical scalar specifying whether path should be overwritten if it already exists.

Value

A directory is produced at path containing the files "matrix.mtx", "barcodes.tsv" and "genes.tsv". A TRUE value is invisibly returned.

Author(s)

Aaron Lun

References

10X Genomics (2017). Gene-Barcode Matrices. https://support.10xgenomics.com/single-cell-gene-expression/software/pipelines/latest/output/matrices

See Also

read10xCounts

Examples

# Mocking up some count data.
library(Matrix)
my.counts <- matrix(rpois(1000, lambda=5), ncol=10, nrow=100)
my.counts <- as(my.counts, "dgCMatrix")
cell.ids <- paste0("BARCODE-", seq_len(ncol(my.counts)))

ngenes <- nrow(my.counts)
gene.ids <- paste0("ENSG0000", seq_len(ngenes))
gene.symb <- paste0("GENE", seq_len(ngenes))

# Writing this to file:
tmpdir <- tempfile()
write10xCounts(tmpdir, my.counts, gene.id=gene.ids, 
    gene.symbol=gene.symb, barcodes=cell.ids)
list.files(tmpdir)

[Package DropletUtils version 1.0.3 Index]