DSBNormalizeProtein | R Documentation |
DSBNormalizeProtein R function: Normalize single cell antibody derived tag (ADT) protein data. This function implements both step I (ambient protein background correction) and step II. (defining and removing cell to cell technical variation) of the dsb normalization method. See <https://www.biorxiv.org/content/10.1101/2020.02.24.963603v3> for details of the algorithm.
DSBNormalizeProtein(
cell_protein_matrix,
empty_drop_matrix,
denoise.counts = TRUE,
use.isotype.control = TRUE,
isotype.control.name.vec = NULL,
define.pseudocount = FALSE,
pseudocount.use,
quantile.clipping = FALSE,
quantile.clip = c(0.001, 0.9995),
scale.factor = c("standardize", "mean.subtract")[1],
return.stats = FALSE
)
cell_protein_matrix |
Raw protein ADT UMI count data to be normalized. Cells - columns Proteins (ADTs) - rows. |
empty_drop_matrix |
Raw empty droplet / background ADT UMI count data used for background correction with Cells - columns and Proteins (ADTs) - rows. This can easily be defined from the raw_feature_bc_matrix output from Cell Ranger or other alignment tools such as kallisto and Cite-Seq-Count. See vignette. |
denoise.counts |
Recommended function default 'denoise.counts = TRUE' and 'use.isotype.control = TRUE'. This runs step II of the dsb algorithm to define and remove cell to cell technical noise. |
use.isotype.control |
Recommended function default 'denoise.counts = TRUE' and 'use.isotype.control = TRUE'. This includes isotype controls in defining the dsb technical component. |
isotype.control.name.vec |
A vector of the names of the isotype control proteins in the rows of the cells and background matrix e.g. isotype.control.name.vec = c('isotype1', 'isotype2'). |
define.pseudocount |
FALSE (default) uses the value 10 optimized for protein ADT data. |
pseudocount.use |
Must be defined if 'define.pseudocount = TRUE'. This is the pseudocount to be added to raw ADT UMI counts. Otherwise the default pseudocount used. |
quantile.clipping |
FALSE (default), if outliers or a large range of values for some proteins are observed (e.g. -50 to 50) these are often from rare outlier cells. re-running the function with 'quantile.clipping = TRUE' will adjust by applying 0.001 and 0.998th quantile value clipping to trim values to those max and min values. If range of normalized values are still very broad and high (e.g. above 40) try setting 'scale.factor = mean.subtract'. |
quantile.clip |
if 'quantile.clipping = TRUE', a vector of the lowest and highest quantiles to clip. These can be tuned to the dataset size. The default c(0.001, 0.9995) optimized to clip only a few of the most extreme outliers. |
scale.factor |
one of 'standardize' or 'mean.subtract'. The recommended default 'standardize' subtracts from the cells the the background droplet matrix mean and divides by the background matrix standard deviation. Values for each protein with this method are interpretable as the number of standard deviations from the mean of the protein background distribution. If 'mean.subtract', subtract the mean without dividing by the standard deviation; can be useful if low background levels detected. |
return.stats |
if TRUE, returns a list, element 1 $dsb_normalized_matrix is the normalized adt matrix element 2 $dsb_stats is the internal stats used by dsb during denoising (the background mean, isotype control values, and the final dsb technical component that is regressed out of the counts) |
Normalized ADT data are returned as a standard R "matrix" of cells (columns), proteins (rows) that can be added to Seurat, SingleCellExperiment or python anndata object - see vignette. If return.stats = TRUE, function returns a list: x$dsb_normalized_matrix normalized matrix, x$protein_stats are mean and sd of log transformed cell, background and the dsb normalized values (as list). x$technical_stats includes the dsb technical component value for each cell and each variable used to calculate the technical component.
Matthew P. Mulè, mattmule@gmail.com
https://doi.org/10.1101/2020.02.24.963603
library(dsb) # load example data cells_citeseq_mtx and empty_drop_matrix included in package
# use a subset of cells and background droplets from example data
cells_citeseq_mtx = cells_citeseq_mtx[ ,1:400]
empty_drop_matrix = empty_drop_citeseq_mtx[ ,1:400]
# example I
adt_norm = dsb::DSBNormalizeProtein(
# step I: remove ambient protein noise reflected in counts from empty droplets
cell_protein_matrix = cells_citeseq_mtx,
empty_drop_matrix = empty_drop_matrix,
# recommended step II: model and remove the technical component of each cell's protein data
denoise.counts = TRUE,
use.isotype.control = TRUE,
isotype.control.name.vec = rownames(cells_citeseq_mtx)[67:70]
)
# example II - experiments without isotype controls
adt_norm = dsb::DSBNormalizeProtein(
cell_protein_matrix = cells_citeseq_mtx,
empty_drop_matrix = empty_drop_matrix,
denoise.counts = FALSE
)
# example III - return dsb internal stats used during denoising for each cell
# returns a 2 element list - the normalized matrix and the internal stats
dsb_output = dsb::DSBNormalizeProtein(
cell_protein_matrix = cells_citeseq_mtx,
empty_drop_matrix = empty_drop_matrix,
isotype.control.name.vec = rownames(cells_citeseq_mtx)[67:70],
return.stats = TRUE
)
# the dsb normalized matrix to be used in downstream analysis is dsb_output$dsb_normalized_matrix
# protein level stats are in dsb_output$protein_stats
# cell-level stats are in dsb_output$technical_stats
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.