read10X: Read 10X Genomics Files

Description Usage Arguments Details Value Author(s) See Also Examples

View source: R/read10X.R

Description

Reads 10X Genomics files containing single-cell RNA-seq UMI counts in Matrix Market format.

Usage

1
read10X(mtx = NULL, genes = NULL, barcodes = NULL, path = ".", DGEList = TRUE)

Arguments

mtx

name of mtx file containing counts in Matrix Exchange Format. Defaults to matrix.mtx or matrix.mtx.gz.

genes

name of file containing gene IDs and names. Defaults to features.tsv or genes.tsv or gzipped versions of the same.

barcodes

optional name of file containing barcodes. Defaults to "barcodes.tsv" or barcodes.tsv.gz.

path

character string giving the directory containing the files. Defaults to the current working directory.

DGEList

logical. If TRUE, a DGEList will be returned, otherwise an unclassed list is returned.

Details

This function reads output files created by the 10X Genomics Cellranger pipeline, see https://support.10xgenomics.com/single-cell-gene-expression/software/pipelines/latest/output/matrices. The UMI counts are assembled into an integer matrix in R with accompanying gene IDs and gene symbols. The results are returned as either a DGEList or an ordinary list.

The files mtx, genes and barcodes can be provided in either gzipped or unzipped versions.

This function creates an ordinary matrix of counts. To read the counts instead into a sparse matrix format, the read10xResults function in the scater package is an alternative.

Value

Either a DGEList object (if DGEList=TRUE) or an ordinary list with the following components:

counts

matrix of counts.

genes

data.frame counting gene symbols.

samples

data.frame containing information about each cell. This will be omitted if barcodes=NULL and DGEList=FALSE.

The only difference between the DGEList or list formats is that the DGEList adds some extra columns to the samples data.frame.

Author(s)

Gordon Smyth

See Also

read10xResults in the scater package.

Examples

1
2
3
4
5
6
7
8
9
## Not run: 
GEO <- "ftp://ftp.ncbi.nlm.nih.gov/geo/samples/GSM2510nnn/GSM2510617/suppl/"
GEOmtx <- paste0(GEO,"GSM2510617_P7-matrix.mtx.gz")
GEOgenes <- paste0(GEO,"GSM2510617_P7-genes.tsv.gz")
download.file(GEOmtx,"matrix.mtx.gz")
download.file(GEOgenes,"genes.tsv.gz")
y <- read10X("matrix.mtx.gz", "genes.tsv.gz")

## End(Not run)

edgeR documentation built on Jan. 16, 2021, 2:03 a.m.