Package 'lfc'

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

Help Index


Subtract the median of the given vector (for normalizing log2 fold changes).

Description

Subtract the median of the given vector (for normalizing log2 fold changes).

Usage

CenterMedian(l)

Arguments

l

Vector of effect sizes

Value

A vector of length 2 containing the two parameters

See Also

PsiLFC

Examples

CenterMedian(rnorm(1000,200))

The log2 fold change distribution

Description

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).

Usage

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)

Arguments

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

Value

The density

Functions

  • dlfc(): Density function

  • plfc(): Distribution function

  • qlfc(): Quantile function

  • rlfc(): random generation

Examples

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

Description

Computes the prior parameters (i.e. pseudocounts incremented by 1) for the log2 fold change distribution

Usage

EmpiricalBayesPrior(A, B, min.sd = 0)

Arguments

A

Vector of counts from condition A

B

Vector of counts from condition B

min.sd

minimal standard deviation of the prior

Value

A vector of length 2 containing the two parameters

See Also

PsiLFC

Examples

EmpiricalBayesPrior(rnorm(1000,200),rnorm(1000,100))

Inverse logit transformation to obtain proportion representation from the log fold change representation.

Description

Inverse logit transformation to obtain proportion representation from the log fold change representation.

Usage

ltop(l)

Arguments

l

Effect size in log2 fold change representation

Value

The proportion representation of the effect size

See Also

ptol

Other Effect size representations: ptol()

Examples

ptol(0)
ptol(1)

Standard LFC effect size estimator

Description

Computes the standard, normalized log2 fold change with given pseudocounts

Usage

NormLFC(A, B, pseudo = c(1, 1), normalizeFun = CenterMedian)

Arguments

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

Value

The vector containing the estimates

Examples

NormLFC(rnorm(1000,200),rnorm(1000,100))

Psi LFC effect size estimator

Description

Computes the optimal effect size estimate and credible intervals if needed.

Usage

PsiLFC(
  A,
  B,
  prior = EmpiricalBayesPrior(A, B),
  normalizeFun = CenterMedian,
  cre = FALSE,
  verbose = FALSE
)

Arguments

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?

Value

Either a vector containing the estimates, or a data frame containing the credible interval as well

Examples

PsiLFC(rnorm(1000,200),rnorm(1000,100))

Psi LFC effect size estimator

Description

Computes the optimal effect size estimate and credible intervals if needed for a Bioconductor SummarizedExperiment object

Usage

PsiLFC.se(se, contrast, cre = FALSE)

Arguments

se

SummarizedExperiment object

contrast

Vector of length 3 (<name>,<A>,<B>)

cre

Compute credible intervals as well? (can also be a vector of quantiles)

Value

Either a vector containing the estimates, or a data frame containing the credible interval as well

Examples

## 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.

Description

Logit transformation to obtain the log fold change representation from the proportion representation.

Usage

ptol(p)

Arguments

p

Effect size in proportion representation

Value

The log2 fold change representation of the effect size

See Also

ltop

Other Effect size representations: ltop()

Examples

ptol(0.5)
ptol(2/3)

Psi LFC effect size estimator for DESeq2

Description

Drop-in replacement for DESeq2's results function for simple settings involving a single variable. Appends the PsiLFC estimate.

Usage

results(object, contrast, cre = FALSE, ...)

Arguments

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

Value

Either a vector containing the estimates, or a data frame containing the credible interval as well