Title: | Log Fold Change Distribution Tools for Working with Ratios of Counts |
---|---|
Description: | Ratios of count data such as obtained from RNA-seq are modelled using Bayesian statistics to derive posteriors for effects sizes. This approach is described in Erhard & Zimmer (2015) <doi:10.1093/nar/gkv696> and Erhard (2018) <doi:10.1093/bioinformatics/bty471>. |
Authors: | Florian Erhard |
Maintainer: | Florian Erhard <[email protected]> |
License: | Apache License (>= 2) |
Version: | 0.2.3 |
Built: | 2025-02-09 05:31:37 UTC |
Source: | https://github.com/erhard-lab/lfc |
Subtract the median of the given vector (for normalizing log2 fold changes).
CenterMedian(l)
CenterMedian(l)
l |
Vector of effect sizes |
A vector of length 2 containing the two parameters
PsiLFC
CenterMedian(rnorm(1000,200))
CenterMedian(rnorm(1000,200))
Density, distribution function, quantile function and random generation for the log2 fold change distribution with parameters ‘a’ and ‘b’ (corresponding to (pseudo-)counts incremented by 1).
dlfc(x, a, b, log = FALSE) plfc(q, a, b, lower.tail = TRUE, log.p = FALSE) qlfc(p, a, b, lower.tail = TRUE, log.p = FALSE) rlfc(n, a, b)
dlfc(x, a, b, log = FALSE) plfc(q, a, b, lower.tail = TRUE, log.p = FALSE) qlfc(p, a, b, lower.tail = TRUE, log.p = FALSE) rlfc(n, a, b)
x , q
|
vector of quantiles |
a |
non-negative parameter |
b |
non-negative parameter |
log , log.p
|
if TRUE, probabilities p are given as log(p) |
lower.tail |
if TRUE (default), probabilities are P[X <= x], otherwise, P[X > x]. |
p |
vector of probabilities |
n |
number of observations |
The density
dlfc()
: Density function
plfc()
: Distribution function
qlfc()
: Quantile function
rlfc()
: random generation
x <- seq (-5,5,by=0.01) plot(x,dlfc(x,1,1))
x <- seq (-5,5,by=0.01) plot(x,dlfc(x,1,1))
Computes the prior parameters (i.e. pseudocounts incremented by 1) for the log2 fold change distribution
EmpiricalBayesPrior(A, B, min.sd = 0)
EmpiricalBayesPrior(A, B, min.sd = 0)
A |
Vector of counts from condition A |
B |
Vector of counts from condition B |
min.sd |
minimal standard deviation of the prior |
A vector of length 2 containing the two parameters
PsiLFC
EmpiricalBayesPrior(rnorm(1000,200),rnorm(1000,100))
EmpiricalBayesPrior(rnorm(1000,200),rnorm(1000,100))
Inverse logit transformation to obtain proportion representation from the log fold change representation.
ltop(l)
ltop(l)
l |
Effect size in log2 fold change representation |
The proportion representation of the effect size
ptol
Other Effect size representations:
ptol()
ptol(0) ptol(1)
ptol(0) ptol(1)
Computes the standard, normalized log2 fold change with given pseudocounts
NormLFC(A, B, pseudo = c(1, 1), normalizeFun = CenterMedian)
NormLFC(A, B, pseudo = c(1, 1), normalizeFun = CenterMedian)
A |
Vector of counts from condition A |
B |
Vector of counts from condition B |
pseudo |
Vector of length 2 of the pseudo counts |
normalizeFun |
Function to normalize the obtained effect sizes |
The vector containing the estimates
NormLFC(rnorm(1000,200),rnorm(1000,100))
NormLFC(rnorm(1000,200),rnorm(1000,100))
Computes the optimal effect size estimate and credible intervals if needed.
PsiLFC( A, B, prior = EmpiricalBayesPrior(A, B), normalizeFun = CenterMedian, cre = FALSE, verbose = FALSE )
PsiLFC( A, B, prior = EmpiricalBayesPrior(A, B), normalizeFun = CenterMedian, cre = FALSE, verbose = FALSE )
A |
Vector of counts from condition A |
B |
Vector of counts from condition B |
prior |
Vector of length 2 of the prior parameters |
normalizeFun |
Function to normalize the obtained effect sizes |
cre |
Compute credible intervals as well? (can also be a vector of quantiles) |
verbose |
verbose status updates? |
Either a vector containing the estimates, or a data frame containing the credible interval as well
PsiLFC(rnorm(1000,200),rnorm(1000,100))
PsiLFC(rnorm(1000,200),rnorm(1000,100))
Computes the optimal effect size estimate and credible intervals if needed for a Bioconductor SummarizedExperiment object
PsiLFC.se(se, contrast, cre = FALSE)
PsiLFC.se(se, contrast, cre = FALSE)
se |
SummarizedExperiment object |
contrast |
Vector of length 3 (<name>,<A>,<B>) |
cre |
Compute credible intervals as well? (can also be a vector of quantiles) |
Either a vector containing the estimates, or a data frame containing the credible interval as well
## Not run: data(airway, package="airway") head(PsiLFC.se(airway,contrast=c("dex","untrt","trt"))) ## End(Not run)
## Not run: data(airway, package="airway") head(PsiLFC.se(airway,contrast=c("dex","untrt","trt"))) ## End(Not run)
Logit transformation to obtain the log fold change representation from the proportion representation.
ptol(p)
ptol(p)
p |
Effect size in proportion representation |
The log2 fold change representation of the effect size
ltop
Other Effect size representations:
ltop()
ptol(0.5) ptol(2/3)
ptol(0.5) ptol(2/3)
Drop-in replacement for DESeq2's results function for simple settings involving a single variable. Appends the PsiLFC estimate.
results(object, contrast, cre = FALSE, ...)
results(object, contrast, cre = FALSE, ...)
object |
the DESeq2DataSet object |
contrast |
Vector of length 3, specifying the variable and the two levels to compute effect sizes for (<name>,<A>,<B>) |
cre |
Compute credible intervals as well? (can also be a vector of quantiles) |
... |
Handed over to DESeq2's results function |
Either a vector containing the estimates, or a data frame containing the credible interval as well