lfcShrink: Shrink log2 fold changes

View source: R/lfcShrink.R

lfcShrinkR Documentation

Shrink log2 fold changes

Description

Adds shrunken log2 fold changes (LFC) and SE to a results table from DESeq run without LFC shrinkage. For consistency with results, the column name lfcSE is used here although what is returned is a posterior SD. Three shrinkage estimators for LFC are available via type (see the vignette for more details on the estimators). The apeglm publication demonstrates that 'apeglm' and 'ashr' outperform the original 'normal' shrinkage estimator.

Usage

lfcShrink(
  dds,
  coef,
  contrast,
  res,
  type = c("apeglm", "ashr", "normal"),
  lfcThreshold = 0,
  svalue = FALSE,
  returnList = FALSE,
  format = c("DataFrame", "GRanges", "GRangesList"),
  saveCols = NULL,
  apeAdapt = TRUE,
  apeMethod = "nbinomCR",
  parallel = FALSE,
  BPPARAM = bpparam(),
  quiet = FALSE,
  ...
)

Arguments

dds

a DESeqDataSet object, after running DESeq

coef

the name or number of the coefficient (LFC) to shrink, consult resultsNames(dds) after running DESeq(dds). note: only coef or contrast can be specified, not both. apeglm requires use of coef. For normal, one of coef or contrast must be provided.

contrast

see argument description in results. only coef or contrast can be specified, not both.

res

a DESeqResults object. Results table produced by the default pipeline, i.e. DESeq followed by results. If not provided, it will be generated internally using coef or contrast. For ashr, if res is provided, then coef and contrast are ignored.

type

"apeglm" is the adaptive Student's t prior shrinkage estimator from the 'apeglm' package; "ashr" is the adaptive shrinkage estimator from the 'ashr' package, using a fitted mixture of normals prior - see the Stephens (2016) reference below for citation; "normal" is the 2014 DESeq2 shrinkage estimator using a Normal prior;

lfcThreshold

a non-negative value which specifies a log2 fold change threshold (as in results). This can be used with any shrinkage type. It will produce new p-values or s-values testing whether the LFC is greater in absolute value than the threshold. The s-values returned in combination with apeglm or ashr provide the probability of FSOS events, "false sign or small", among the tests with equal or smaller s-value than a given gene's s-value, where "small" is specified by lfcThreshold.

svalue

logical, should p-values and adjusted p-values be replaced with s-values when using apeglm or ashr. s-values provide the probability of false signs among the tests with equal or smaller s-value than a given given's s-value. See Stephens (2016) reference on s-values.

returnList

logical, should lfcShrink return a list, where the first element is the results table, and the second element is the output of apeglm or ashr

format

same as defined in results, either "DataFrame", "GRanges", or "GRangesList"

saveCols

same as defined in results, which metadata columns to pass into output

apeAdapt

logical, should apeglm use the MLE estimates of LFC to adapt the prior, or use default or specified prior.control

apeMethod

what method to run apeglm, which can differ in terms of speed

parallel

if FALSE, no parallelization. if TRUE, parallel execution using BiocParallel, see same argument of DESeq parallelization only used with normal or apeglm

BPPARAM

see same argument of DESeq

quiet

whether to print messages

...

arguments passed to apeglm and ashr

Details

See vignette for a comparison of shrinkage estimators on an example dataset. For all shrinkage methods, details on the prior is included in priorInfo(res), including the fitted_g mixture for ashr.

For type="apeglm": Specifying apeglm passes along DESeq2 MLE log2 fold changes and standard errors to the apeglm function in the apeglm package, and re-estimates posterior LFCs for the coefficient specified by coef.

For type="ashr": Specifying ashr passes along DESeq2 MLE log2 fold changes and standard errors to the ash function in the ashr package, with arguments mixcompdist="normal" and method="shrink".

For type="normal": For design as a formula, shrinkage cannot be applied to coefficients in a model with interaction terms. For user-supplied model matrices, shrinkage is only supported via coef. coef will make use of standard model matrices, while contrast will make use of expanded model matrices, and for the latter, a results object should be provided to res. With numeric- or list-style contrasts, it is possible to use lfcShrink, but likely easier to use DESeq with betaPrior=TRUE followed by results, because the numeric or list should reference the coefficients from the expanded model matrix. These coefficients will be printed to console if 'contrast' is used.

Value

a DESeqResults object with the log2FoldChange and lfcSE columns replaced with shrunken LFC and SE. For consistency with results (and similar to the output of bayesglm) the column name lfcSE is used here, although what is returned is a posterior SD. For normal and for apeglm the estimate is the posterior mode, for ashr it is the posterior mean. priorInfo(res) contains information about the shrinkage procedure, relevant to the various methods specified by type.

References

Publications for the following shrinkage estimators:

type="apeglm":

Zhu, A., Ibrahim, J.G., Love, M.I. (2018) Heavy-tailed prior distributions for sequence count data: removing the noise and preserving large differences. Bioinformatics. https://doi.org/10.1093/bioinformatics/bty895

type="ashr":

Stephens, M. (2016) False discovery rates: a new deal. Biostatistics, 18:2. https://doi.org/10.1093/biostatistics/kxw041

type="normal":

Love, M.I., Huber, W., Anders, S. (2014) Moderated estimation of fold change and dispersion for RNA-seq data with DESeq2. Genome Biology, 15:550. https://doi.org/10.1186/s13059-014-0550-8

Related work, the bayesglm function in the arm package:

Gelman, A., Jakulin, A., Pittau, M.G. and Su, Y.-S. (2009). A Weakly Informative Default Prior Distribution For Logistic And Other Regression Models. The Annals of Applied Statistics 2:4. http://www.stat.columbia.edu/~gelman/research/published/ priors11.pdf

Examples


set.seed(1)
dds <- makeExampleDESeqDataSet(n=500,betaSD=1)
dds <- DESeq(dds)
res <- results(dds)

# these are the coefficients from the model
# we can specify them using 'coef' by name or number below
resultsNames(dds)

res.ape <- lfcShrink(dds=dds, coef=2, type="apeglm")
res.ash <- lfcShrink(dds=dds, coef=2, type="ashr")
res.norm <- lfcShrink(dds=dds, coef=2, type="normal")


mikelove/DESeq2 documentation built on July 25, 2024, 11:11 p.m.