Last data update: 2014.03.03

R: Normalization by deconvolution
Deconvolution MethodsR Documentation

Normalization by deconvolution

Description

Methods to normalize single-cell RNA-seq data by deconvolving size factors from cell pools.

Usage

## S4 method for signature 'matrix'
computeSumFactors(x, sizes=c(20, 40, 60, 80, 100), 
    clusters=NULL, ref.clust=NULL, positive=FALSE, errors=FALSE)
## S4 method for signature 'SCESet'
computeSumFactors(x, ..., assay="counts", get.spikes=FALSE)

Arguments

x

A numeric count matrix where rows are genes and columns are cells. Alternatively, a SCESet object containing such a matrix.

sizes

A numeric vector of pool sizes, i.e., number of cells per pool.

clusters

An optional factor specifying which cells belong to which cluster, for deconvolution within clusters.

ref.clust

A level of clusters to be used as the reference cluster for inter-cluster normalization.

positive

A logical scalar indicating whether linear inverse models should be used to enforce positive estimates.

errors

A logical scalar indicating whether the standard error should be returned.

...

Additional arguments to pass to computeSumFactors,matrix-method.

assay

A string specifying which assay values to use, e.g., counts or exprs.

get.spikes

A logical scalar specifying whether spike-in transcripts should be used.

Value

For computeSumFactors,matrix-method, a numeric vector of size factors for all cells in x is returned.

For computeSumFactors,SCESet-method, an object of class x is returned containing the vector of size factors in sizeFactors(x).

If errors=TRUE, the standard error of the size factor estimates is stored as the "standard.error" field of the attributes of the returned vector.

Overview of the deconvolution method

The computeSumFactors function provides an implementation of the deconvolution strategy for normalization. Briefly, a pool of cells is selected and the counts for those cells are summed together. The count sums for this pool is normalized against an average reference pseudo-cell, constructed by averaging the counts across all cells. This defines a size factor for the pool as the median ratio between the count sums and the average across all genes.

Now, the bias for the pool is equal to the sum of the biases for the constituent cells. The same applies for the size factors (which are effectively estimates of the bias for each cell). This means that the size factor for the pool can be written as a linear equation of the size factors for the cells. Repeating this process for multiple pools will yield a linear system that can be solved to obtain the size factors for the individual cells.

In this manner, pool-based factors are deconvolved to yield the relevant cell-based factors. The advantage is that the pool-based estimates are more accurate, as summation reduces the number of stochastic zeroes and the associated bias of the size factor estimate. This accuracy will feed back into the deconvolution process, thus improving the accuracy of the cell-based size factors. The standard error of the estimates can be obtained by setting errors=TRUE.

Normalization within and between clusters

In general, it is more appropriate to pool more similar cells to avoid violating the assumption of a non-DE majority of genes across the data set. This can be done by specifying the clusters argument where cells in each cluster have similar expression profiles. Deconvolution is subsequently applied on the cells within each cluster. Each cluster should contain a sufficient number of cells for pooling – twice the maximum value of sizes is recommended. A convenince function quickCluster is provided for rapid clustering based on Spearman's rank correlation.

Size factors computed within each cluster must be rescaled for comparison between clusters. This is done by normalizing between clusters to identify the rescaling factor. One cluster is chosen as a “reference” (by default, that with the median of the mean per-cell library sizes is used) to which all others are normalized. Ideally, a cluster that is not extremely different from all other clusters should be used as the reference. This can be specified using ref.clust if there is prior knowledge about which cluster is most suitable, e.g., from PCA or t-SNE plots.

Additional details about pooling and deconvolution

Within each cluster (if not specified, all cells are put into a single cluster), cells are sorted by increasing library size and a sliding window is applied to this ordering. Each location of the window defines a cell pool with similar library sizes. This avoids inflated estimation errors for very small cells when they are pooled with very large cells. Sliding the window will construct a linear system. This is repeated with all window sizes in sizes to obtain an over-determined system that can be solved with methods like the QR decomposition.

In theory, it is possible to obtain negative estimates for the size factors. These are most likely for very small library sizes and are obviously nonsensical. Some protection can be provided by setting positive=TRUE, which will use linear inverse models to solve the system. This ensures that non-negative values for the size factors will always be obtained. Note that some cells may still have size factors of zero and should be removed prior to downstream analysis. Such occurrences are unavoidable – rather, the aim is to prevent negative values from affecting the estimates for all other cells.

By default, get.spikes=FALSE which means that spike-in transcripts are not included in the set of genes used for deconvolution. This is because they can behave differently from the endogenous genes. Users wanting to perform spike-in normalization should see computeSpikeFactors instead.

Author(s)

Aaron Lun and Karsten Bach

See Also

quickCluster

Examples

set.seed(100)
popsize <- 800
ngenes <- 10000
all.facs <- 2^rnorm(popsize, sd=0.5)
counts <- matrix(rnbinom(ngenes*popsize, mu=all.facs*10, size=1), ncol=popsize, byrow=TRUE)
out.facs <- computeSumFactors(counts)

Results


R version 3.3.1 (2016-06-21) -- "Bug in Your Hair"
Copyright (C) 2016 The R Foundation for Statistical Computing
Platform: x86_64-pc-linux-gnu (64-bit)

R is free software and comes with ABSOLUTELY NO WARRANTY.
You are welcome to redistribute it under certain conditions.
Type 'license()' or 'licence()' for distribution details.

R is a collaborative project with many contributors.
Type 'contributors()' for more information and
'citation()' on how to cite R or R packages in publications.

Type 'demo()' for some demos, 'help()' for on-line help, or
'help.start()' for an HTML browser interface to help.
Type 'q()' to quit R.

> library(scran)
Loading required package: BiocParallel
Loading required package: scater
Loading required package: Biobase
Loading required package: BiocGenerics
Loading required package: parallel

Attaching package: 'BiocGenerics'

The following objects are masked from 'package:parallel':

    clusterApply, clusterApplyLB, clusterCall, clusterEvalQ,
    clusterExport, clusterMap, parApply, parCapply, parLapply,
    parLapplyLB, parRapply, parSapply, parSapplyLB

The following objects are masked from 'package:stats':

    IQR, mad, xtabs

The following objects are masked from 'package:base':

    Filter, Find, Map, Position, Reduce, anyDuplicated, append,
    as.data.frame, cbind, colnames, do.call, duplicated, eval, evalq,
    get, grep, grepl, intersect, is.unsorted, lapply, lengths, mapply,
    match, mget, order, paste, pmax, pmax.int, pmin, pmin.int, rank,
    rbind, rownames, sapply, setdiff, sort, table, tapply, union,
    unique, unsplit

Welcome to Bioconductor

    Vignettes contain introductory material; view with
    'browseVignettes()'. To cite Bioconductor, see
    'citation("Biobase")', and for packages 'citation("pkgname")'.

Loading required package: ggplot2

Attaching package: 'scater'

The following object is masked from 'package:stats':

    filter

> png(filename="/home/ddbj/snapshot/RGM3/R_BC/result/scran/computeSumFactors.Rd_%03d_medium.png", width=480, height=480)
> ### Name: Deconvolution Methods
> ### Title: Normalization by deconvolution
> ### Aliases: computeSumFactors computeSumFactors,matrix-method
> ###   computeSumFactors,SCESet-method
> ### Keywords: normalization
> 
> ### ** Examples
> 
> set.seed(100)
> popsize <- 800
> ngenes <- 10000
> all.facs <- 2^rnorm(popsize, sd=0.5)
> counts <- matrix(rnbinom(ngenes*popsize, mu=all.facs*10, size=1), ncol=popsize, byrow=TRUE)
> out.facs <- computeSumFactors(counts)
> 
> 
> 
> 
> 
> dev.off()
null device 
          1 
>