lfcShrink {DESeq2}R Documentation

Shrink log2 fold changes

Description

Adds shrunken log2 fold changes (LFC) and SE to a results table from DESeq run without LFC shrinkage. Three shrinkage esimators for LFC are available via type.

Usage

lfcShrink(dds, coef, contrast, res, type = c("normal", "apeglm", "ashr"),
  lfcThreshold = 0, svalue = FALSE, returnList = FALSE, 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.

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

"normal" is the original DESeq2 shrinkage estimator; "apeglm" is the adaptive 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

lfcThreshold

a non-negative value which specifies a log2 fold change threshold (as in results). This can be used in conjunction with normal and apeglm, where 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 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

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

As of DESeq2 version 1.18, type="apeglm" and type="ashr" are new features, and still under development. 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. 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". 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. The integration of shrinkage methods from external packages will likely evolve over time. We will likely incorporate an lfcThreshold argument which can be passed to apeglm to specify regions of the posterior at an arbitrary threshold.

For normal, and design as a formula, shrinkage cannot be applied to coefficients in a model with interaction terms. For normal and user-supplied model matrices, shrinkage is only supported via coef. For normal 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 with normal.

Value

a DESeqResults object with the log2FoldChange and lfcSE columns replaced with shrunken LFC and SE. priorInfo(res) contains information about the shrinkage procedure, relevant to the various methods specified by type.

References

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

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. bioRxiv. https://doi.org/10.1101/303255

type="ashr":

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

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.shr <- lfcShrink(dds=dds, coef=2)
res.shr <- lfcShrink(dds=dds, contrast=c("condition","B","A"))
res.ape <- lfcShrink(dds=dds, coef=2, type="apeglm")
res.ash <- lfcShrink(dds=dds, coef=2, type="ashr")


[Package DESeq2 version 1.20.0 Index]