Skip to contents

Wrapper function that helps set up DESeq2::lfcShrink() to shrink LFC values for a pairwise contrast via apeglm, without having to manually relevel factor reference levels to use the required coef argument.

Usage

apeglmResults(object, ...)

# S4 method for DESeqDataSet
apeglmResults(object, contrast, res, ...)

Arguments

object

Object.

...

Additional arguments.

contrast

character(3). Pairwise contrast vector:

  1. factor: Grouping factor. Corresponds to column name in colData().

  2. numerator: Numerator samples.

  3. denominator: Denominator samples.

Numerator and denominator values correspond to grouping factor column. See results() for details. Note that we're intentionally being more strict about the input format here.

res

DESeqResults. Results containing unshrunken LFC values, generated with results().

Value

DESeqResults, with apeglm adaptive shrinkage applied to fold change values.

Details

Dynamically sets reference factor levels, as recommended by DESeq2 vignette. Matches contrast input internally to corresponding coef corresponding to values in resultsNames().

Runs nbinomWaldTest() via DESeq(), followed by lfcShrink().

Note

Updated 2020-08-17.

See also

Examples

## DESeqDataSet ====
if (requireNamespace("apeglm", quietly = TRUE)) {
    dds <- DESeq2::makeExampleDESeqDataSet(n = 1000L, m = 12L)
    dds$condition <- factor(rep(LETTERS[seq_len(4L)], each = 3L))
    dds <- DESeq2::DESeq(dds)
    resultsNames(dds)

    ## Contrast C vs. B.
    contrast <- c(factor = "condition", numerator = "C", denominator = "B")

    ## Unshrunken DESeqResults.
    res <- DESeq2::results(dds, contrast = contrast)
    class(res)
    lfcShrinkType(res)

    ## Shrunken DESeqResults, using apeglm via `lfcShrink()`.
    shrink <- apeglmResults(
        object = dds,
        contrast = contrast,
        res = res
    )
    class(shrink)
    lfcShrinkType(shrink)
}
#> estimating size factors
#> estimating dispersions
#> gene-wise dispersion estimates
#> mean-dispersion relationship
#> final dispersion estimates
#> fitting model and testing
#> design: ~condition
#> using pre-existing size factors
#> estimating dispersions
#> gene-wise dispersion estimates: 6 workers
#> mean-dispersion relationship
#> final dispersion estimates, fitting model and testing: 6 workers
#> contrast: condition_C_vs_B
#> coef: 3
#> [1] "apeglm"