View source: R/delta_method_transform.R
acosh_transform | R Documentation |
Delta method-based variance stabilizing transformation
acosh_transform(
data,
overdispersion = 0.05,
size_factors = TRUE,
...,
on_disk = NULL,
verbose = FALSE
)
shifted_log_transform(
data,
overdispersion = 0.05,
pseudo_count = 1/(4 * overdispersion),
size_factors = TRUE,
minimum_overdispersion = 0.001,
...,
on_disk = NULL,
verbose = FALSE
)
data |
any matrix-like object (e.g. matrix, dgCMatrix, DelayedArray, HDF5Matrix)
with one column per sample and row per gene. It can also be an object of type |
overdispersion |
the simplest count model is the Poisson model. However, the Poisson model
assumes that
Note that |
size_factors |
in large scale experiments, each sample is typically of different size
(for example different sequencing depths). A size factor is an internal mechanism of GLMs to
correct for this effect. |
... |
additional parameters for |
on_disk |
a boolean that indicates if the dataset is loaded into memory or if it is kept on disk
to reduce the memory usage. Processing in memory can be significantly faster than on disk.
Default: |
verbose |
boolean that decides if information about the individual steps are printed.
Default: |
pseudo_count |
instead of specifying the overdispersion, the
|
minimum_overdispersion |
the |
a matrix (or a vector if the input is a vector) with the transformed values.
acosh_transform
: 1/sqrt(alpha)
acosh(2 * alpha * x + 1)
shifted_log_transform
: 1/sqrt(alpha) log(4 * alpha * x + 1)
Ahlmann-Eltze, Constantin, and Wolfgang Huber. "Transformation and Preprocessing of Single-Cell RNA-Seq Data." bioRxiv (2021).
Ahlmann-Eltze, Constantin, and Wolfgang Huber. "glmGamPoi: Fitting Gamma-Poisson Generalized Linear Models on Single Cell Count Data." Bioinformatics (2020)
Dunn, Peter K., and Gordon K. Smyth. "Randomized quantile residuals." Journal of Computational and Graphical Statistics 5.3 (1996): 236-244.
Hafemeister, Christoph, and Rahul Satija. "Normalization and variance stabilization of single-cell RNA-seq data using regularized negative binomial regression." Genome biology 20.1 (2019): 1-15.
Hafemeister, Christoph, and Rahul Satija. "Analyzing scRNA-seq data with the sctransform and offset models." (2020)
Lause, Jan, Philipp Berens, and Dmitry Kobak. "Analytic Pearson residuals for normalization of single-cell RNA-seq UMI data." Genome Biology (2021).
acosh_transform
, shifted_log_transform
, and residual_transform
# Load a single cell dataset
sce <- TENxPBMCData::TENxPBMCData("pbmc4k")
# Reduce size for this example
set.seed(1)
sce_red <- sce[sample(which(rowSums2(counts(sce)) > 0), 1000),
sample(ncol(sce), 200)]
assay(sce_red, "acosh") <- acosh_transform(sce_red)
assay(sce_red, "shifted_log") <- shifted_log_transform(sce_red)
plot(rowMeans2(assay(sce_red, "acosh")), rowVars(assay(sce_red, "acosh")), log = "x")
points(rowMeans2(assay(sce_red, "shifted_log")), rowVars(assay(sce_red, "shifted_log")),
col = "red")
# Sqrt transformation
sqrt_dat <- acosh_transform(sce_red, overdispersion = 0, size_factor = 1)
plot(2 * sqrt(assay(sce_red))[,1], sqrt_dat[,1]); abline(0,1)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.