nullResiduals | R Documentation |
Computes deviance or Pearson residuals for count data based on a multinomial null model that assumes each feature has a constant rate. The residuals matrix can be analyzed with standard PCA as a fast approximation to GLM-PCA.
nullResiduals(object, ...)
## S4 method for signature 'SummarizedExperiment'
nullResiduals(
object,
assay = "counts",
fam = c("binomial", "poisson"),
type = c("deviance", "pearson"),
batch = NULL
)
## S4 method for signature 'SingleCellExperiment'
nullResiduals(
object,
assay = "counts",
fam = c("binomial", "poisson"),
type = c("deviance", "pearson"),
batch = NULL
)
## S4 method for signature 'matrix'
nullResiduals(
object,
fam = c("binomial", "poisson"),
type = c("deviance", "pearson"),
batch = NULL
)
## S4 method for signature 'Matrix'
nullResiduals(
object,
fam = c("binomial", "poisson"),
type = c("deviance", "pearson"),
batch = NULL
)
## S4 method for signature 'ANY'
nullResiduals(
object,
fam = c("binomial", "poisson"),
type = c("deviance", "pearson"),
batch = NULL
)
object |
The object on which to compute residuals. It can be a
matrix-like object (e.g. matrix, Matrix, DelayedMatrix, HDF5Matrix) with
genes in the rows and samples in the columns. Specialized methods are
defined for objects inheriting from SummarizedExperiment (such as
|
... |
for the generic, additional arguments to pass to object-specific methods. |
assay |
a string or integer specifying which assay contains the count
data (default = 'counts'). Ignored if |
fam |
a string specifying the model type to be used for calculating the residuals. Binomial (the default) is the closest approximation to multinomial, but Poisson may be faster to compute and often is very similar to binomial. |
type |
should deviance or Pearson residuals be used? |
batch |
an optional factor indicating batch membership of observations. If provided, the null model is computed within each batch separately to regress out the batch effect from the resulting residuals. |
This function should be used only on the un-normalized counts. It was originally designed for single-cell RNA-seq counts obtained by the use of unique molecular identifiers (UMIs) and has not been tested on read count data without UMIs or other data types.
Note that even though sparse Matrix objects are accepted as input, they are internally coerced to dense matrix before processing, because the output is always a dense matrix since the residuals transformation is not sparsity preserving. To avoid memory issues, it is recommended to perform feature selection first and subset the number of features to a smaller size prior to computing the residuals.
The original SingleCellExperiment
or
SummarizedExperiment
object with the residuals appended as a new
assay. The assay name will be fam_type_residuals (eg,
binomial_deviance_residuals). If the input was a matrix, output is a dense
matrix containing the residuals.
Townes FW, Hicks SC, Aryee MJ, and Irizarry RA (2019). Feature Selection and Dimension Reduction for Single Cell RNA-Seq based on a Multinomial Model. Genome Biology https://doi.org/10.1186/s13059-019-1861-6
ncells <- 100
u <- matrix(rpois(20000, 5), ncol=ncells)
sce <- SingleCellExperiment::SingleCellExperiment(assays=list(counts=u))
nullResiduals(sce)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.