defaultDrops | R Documentation |
Call cells according to the number of UMIs associated with each barcode, as implemented in CellRanger version 2.
defaultDrops(m, ...)
## S4 method for signature 'ANY'
defaultDrops(m, expected = 3000, upper.quant = 0.99, lower.prop = 0.1)
## S4 method for signature 'SummarizedExperiment'
defaultDrops(m, ..., assay.type = "counts")
m |
A numeric matrix-like object containing counts, where columns represent barcoded droplets and rows represent features. The matrix should only contain barcodes for an individual sample, prior to any filtering for cells. Alternatively, a SummarizedExperiment containing such a matrix. |
... |
For the generic, further arguments to pass to individual methods. For the SummarizedExperiment method, further arguments to pass to the ANY method. |
expected |
A numeric scalar specifying the expected number of cells in this sample, as specified in the call to CellRanger. |
upper.quant |
A numeric scalar between 0 and 1 specifying the quantile of the top |
lower.prop |
A numeric scalar between 0 and 1 specifying the fraction of molecules of the |
assay.type |
Integer or string specifying the assay containing the count matrix. |
The defaultDrops
function will call cells based on library size similarly to the CellRanger software suite from 10X Genomics.
Default arguments correspond to an exact reproduction of CellRanger's algorithm, where the number of expected cells was also left at CellRanger default value.
The method computes the upper.quant
quantile of the top expected
barcodes, ordered by decreasing number of UMIs.
Any barcodes containing more molecules than lower.prop
times this quantile is considered to be a cell, and is retained for further analysis.
This method may be vulnerable to calling very well-captured background RNA as cells, or missing real cells with low RNA content.
See ?emptyDrops
for an alternative approach for cell calling.
A logical vector of length ncol(m)
, indicating whether each column of m
was called as a cell.
Jonathan Griffiths
10X Genomics (2017). Cell Ranger Algorithms Overview. https://support.10xgenomics.com/single-cell-gene-expression/software/pipelines/latest/algorithms/overview
emptyDrops
, for another method for calling cells.
# Mocking up some data:
set.seed(0)
my.counts <- DropletUtils:::simCounts()
# Identify likely cell-containing droplets.
called <- defaultDrops(my.counts)
table(called)
# Get matrix of called cells.
cell.counts <- my.counts[, called]
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.