Package 'iSEEde'

Title: iSEE extension for panels related to differential expression analysis
Description: This package contains diverse functionality to extend the usage of the iSEE package, including additional classes for the panels or modes facilitating the analysis of differential expression results. This package does not perform differential expression. Instead, it provides methods to embed precomputed differential expression results in a SummarizedExperiment object, in a manner that is compatible with interactive visualisation in iSEE applications.
Authors: Kevin Rue-Albrecht [aut, cre] , Thomas Sandmann [ctb] , Denali Therapeutics [fnd]
Maintainer: Kevin Rue-Albrecht <[email protected]>
License: Artistic-2.0
Version: 1.3.1
Built: 2024-10-02 03:27:44 UTC
Source: https://github.com/bioc/iSEEde

Help Index


Extract contrast results embedded in a SummarizedExperiment object

Description

contrastResults returns either all contrasts results stored in object or a single contrast result by name.

contrastResultsNames returns the names of contrast results embedded in object.

Usage

contrastResultsNames(object)

contrastResults(object, name)

Arguments

object

A SummarizedExperiment object.

name

(Optional) Name of a single contrast result name to extract. Use contrastResultsNames(object) to list available names.

Value

For contrastResultsNames: the names of embedded contrast results available.

For contrastResults: a DataFrame of differential expression statistics.

If name is missing, contrastResults returns a nested DataFrame in which each column contains the results of a single contrast. If name is given, contrastResults returns a DataFrame that contains the results of a single contrast.

Examples

library("iSEEde")
library("airway")
library("DESeq2")
library("iSEE")

##
# Example data ----
##

data("airway")
airway$dex <- relevel(airway$dex, "untrt")

dds <- DESeqDataSet(airway, ~ 0 + dex + cell)

dds <- DESeq(dds)
res_deseq2 <- results(dds, contrast = list("dextrt", "dexuntrt"))
airway <- embedContrastResults(res_deseq2, airway, name = "dex: trt vs untrt")

##
# List result names ---
##

contrastResultsNames(airway)

##
# Extract results ---
##

contrastResults(airway)
contrastResults(airway, "dex: trt vs untrt")

Generics for Differential Expression Results

Description

An overview of the generics for accessing common pieces of information in differential expression results.

Definitions

  • pValue(x) returns a named numeric vector of raw p-values.

  • log2FoldChange(x) returns a named numeric vector of log2-fold-change values.

  • averageLog2(x) returns a named numeric vector of average log2-expression values.

Author(s)

Kevin Rue-Albrecht

Examples

showMethods(pValue)
showMethods(log2FoldChange)
showMethods(averageLog2)

The DETable class

Description

The DETable class is a RowTable subclass that is dedicated to creating a volcano plot. It retrieves the table of results for the selected differential expression contrast and creates an interactive table where each row represents a feature.

Slot overview

The following slots control the test procedure:

  • ContrastName, a character scalar indicating the name of the contrast to display.

  • RoundDigits, a logical scalar indicating whether to round numeric values (see SignifDigits).

  • SignifDigits, an integer scalar indicating the number of significant digits to use for rounding numbers (see RoundDigits).

In addition, this class inherits all slots from its parent RowTable and Table classes.

Examples

x <- DETable()
x

The iSEEDESeq2Results class

Description

The iSEEDESeq2Results class is used to provide a common interface to differential expression results produced by the DESeq2 package. It provides methods to access common differential expression statistics (e.g., log2 fold-change, p-value, log2 average abundance).

Details

This class inherits all its slots directly from its parent class DataFrame.

Constructor

iSEEDESeq2Results(data, row.names = rownames(data)) creates an instance of a iSEEDESeq2Results class, with:

data

A data.frame produced by DESeq2::results() or DESeq2::lfcShrink().

row.names

The character vector of rownames for the SummarizedExperiment object in which the object is to be embedded. Must be a superset of rownames(data).

Supported methods

  • embedContrastResults(x, se, name, ...) embeds x in se under the identifier name. See embedContrastResults() for more details.

  • pValue(x) returns the vector of raw p-values.

  • log2FoldChange(x) returns the vector of log2-fold-change values.

  • averageLog2(x) returns the vector of average log2-expression values.

Author(s)

Kevin Rue-Albrecht

Examples

library(DESeq2)

##
# From DESeq2::DESeq() ----
##

cnts <- matrix(rnbinom(n = 1000, mu = 100, size = 1 / 0.5), ncol = 10)
rownames(cnts) <- paste("Gene", 1:100)
cond <- factor(rep(1:2, each = 5))

# object construction
dds <- DESeqDataSetFromMatrix(cnts, DataFrame(cond), ~cond)

# standard analysis
dds <- DESeq(dds)
res <- results(dds)
head(res)

##
# iSEEDESeq2Results ----
##

# Embed the DESeq2 results in the SummarizedExperiment object
dds <- embedContrastResults(res, dds, name = "DESeq2")

##
# Access ----
##

contrastResultsNames(dds)
contrastResults(dds)
contrastResults(dds, "DESeq2")

head(pValue(contrastResults(dds, "DESeq2")))
head(log2FoldChange(contrastResults(dds, "DESeq2")))
head(averageLog2(contrastResults(dds, "DESeq2")))

The iSEEedgeRResults class

Description

The iSEEedgeRResults class is used to provide a common interface to differential expression results produced by the edgeR package. It provides methods to access common differential expression statistics (e.g., log fold-change, p-value, log2 average abundance).

Details

This class inherits all its slots directly from its parent class DataFrame.

Constructor

iSEEedgeRResults(data, row.names = rownames(data)) creates an instance of a iSEEedgeRResults class, with:

data

A data.frame produced by edgeR::topTags().

row.names

The character vector of rownames for the SummarizedExperiment object in which the object is to be embedded. Must be a superset of rownames(data).

Supported methods

  • embedContrastResults(x, se, name, ...) embeds x in se under the identifier name. See embedContrastResults() for more details.

  • pValue(x) returns the vector of raw p-values.

  • log2FoldChange(x) returns the vector of log2-fold-change values.

  • averageLog2(x) returns the vector of average log2-expression values.

Author(s)

Kevin Rue-Albrecht

Examples

library(edgeR)
library(SummarizedExperiment)

##
# From edgeR::glmLRT() ----
##

nlibs <- 3
ngenes <- 100
dispersion.true <- 0.1

# Make first gene respond to covariate x
x <- 0:2
design <- model.matrix(~x)
beta.true <- cbind(Beta1=2,Beta2=c(2,rep(0,ngenes-1)))
mu.true <- 2^(beta.true %*% t(design))

# Generate count data
y <- rnbinom(ngenes*nlibs,mu=mu.true,size=1/dispersion.true)
y <- matrix(y,ngenes,nlibs)
colnames(y) <- c("x0","x1","x2")
rownames(y) <- paste("gene",1:ngenes,sep=".")
d <- DGEList(y)

# Normalize
d <- calcNormFactors(d)

# Fit the NB GLMs
fit <- glmFit(d, design, dispersion=dispersion.true)

# Likelihood ratio tests for trend
results <- glmLRT(fit, coef=2)
tt <- topTags(results)

##
# iSEEedgeRResults ----
##

# Simulate the original SummarizedExperiment object
se <- SummarizedExperiment(assays = list(counts = d$counts))

# Embed the edgeR results in the SummarizedExperiment object
se <- embedContrastResults(tt, se, name = "edgeR")

##
# Access ----
##

contrastResultsNames(se)
contrastResults(se)
contrastResults(se, "edgeR")

head(pValue(contrastResults(se, "edgeR")))
head(log2FoldChange(contrastResults(se, "edgeR")))
head(averageLog2(contrastResults(se, "edgeR")))

The iSEELimmaResults class

Description

The iSEELimmaResults class is used to provide a common interface to differential expression results produced by the limma package. It provides methods to access common differential expression statistics (e.g., log fold-change, p-value, log2 average abundance).

Details

This class inherits all its slots directly from its parent class DataFrame.

Constructor

iSEELimmaResults(data, row.names = rownames(data)) creates an instance of a iSEELimmaResults class, with:

data

A data.frame produced by limma::topTable().

row.names

The character vector of rownames for the SummarizedExperiment object in which the object is to be embedded. Must be a superset of rownames(data).

Supported methods

  • embedContrastResults(x, se, name, class = "limma", ...) embeds x in se under the identifier name. See embedContrastResults() for more details.

  • pValue(x) returns the vector of raw p-values.

  • log2FoldChange(x) returns the vector of log2-fold-change values.

  • averageLog2(x) returns the vector of average log2-expression values.

Author(s)

Kevin Rue-Albrecht

Examples

library(limma)
library(SummarizedExperiment)

##
# From limma::lmFit() ----
##

sd <- 0.3 * sqrt(4 / rchisq(100, df = 4))
y <- matrix(rnorm(100 * 6, sd = sd), 100, 6)
rownames(y) <- paste("Gene", 1:100)
y[1:2, 4:6] <- y[1:2, 4:6] + 2
design <- cbind(Grp1 = 1, Grp2vs1 = c(0, 0, 0, 1, 1, 1))

fit <- lmFit(y, design)
fit <- eBayes(fit)
tt <- topTable(fit, coef = 2)
head(tt)

##
# iSEELimmaResults ----
##

# Simulate the original SummarizedExperiment object
se <- SummarizedExperiment(assays = list(counts = y))

# Embed the Limma-Voom results in the SummarizedExperiment object
se <- embedContrastResults(tt, se, name = "Limma-Voom", class = "limma")

##
# Access ----
##

contrastResultsNames(se)
contrastResults(se)
contrastResults(se, "Limma-Voom")

head(pValue(contrastResults(se, "Limma-Voom")))
head(log2FoldChange(contrastResults(se, "Limma-Voom")))
head(averageLog2(contrastResults(se, "Limma-Voom")))

The LogFCLogFCPlot class

Description

The LogFCLogFCPlot class is a RowDataPlot subclass that is dedicated to comparing the log-fold-change value of two contrasts. It retrieves the log-fold change of the two selected contrasts and creates a row-based plot where each point represents a feature.

Slot overview

The following slots control the test procedure:

  • ContrastNameX, a character scalar indicating the name of the contrast to display on the x-axis.

  • ContrastNameY, a character scalar indicating the name of the contrast to display on the y-axis.

In addition, this class inherits all slots from its parent RowDotPlot, DotPlot, and Panel classes.

Examples

x <- LogFCLogFCPlot()
x

The MAPlot class

Description

The MAPlot is a RowDataPlot subclass that is dedicated to creating an MA plot. It retrieves the log-fold change (M) and mean average (A) values and creates a row-based plot where each point represents a feature.

Slot overview

The following slots control the test procedure:

  • ContrastName, a character scalar indicating the name of the contrast to display.

In addition, this class inherits all slots from its parent RowDotPlot, DotPlot, and Panel classes.

Examples

x <- MAPlot()
x

Generics for Embbedding Results into a SummarizedExperiment Object

Description

An overview of the generics for embedding results into a SummarizedExperiment object, in a format compatible with iSEEde.

Usage

embedContrastResults(x, se, name, ...)

embedContrastResultsMethods

## S4 method for signature 'ANY'
embedContrastResults(x, se, name, ...)

## S4 method for signature 'data.frame'
embedContrastResults(x, se, name, class, ...)

Arguments

x

Object to be embedded.

se

A SummarizedExperiment object.

name

Identifier for the embedded object.

...

Arguments passed to and from other methods.

class

Class to use for embedding x. Only used when class(x) does not uniquely identify the package that generated the object.

Format

embedContrastResultsMethods: Named character vector mapping keywords to class names designed to store differential expression results.

Value

An updated SummarizedExperiment object that contains the embedded object.

Definitions

Author(s)

Kevin Rue-Albrecht

Examples

embedContrastResultsMethods

showMethods(embedContrastResults)

The VolcanoPlot class

Description

The VolcanoPlot is a RowDataPlot subclass that is dedicated to creating a volcano plot. It retrieves the log-fold change and p-value from and creates a row-based plot where each point represents a feature.

Slot overview

The following slots control the test procedure:

  • ContrastName, a character scalar indicating the name of the contrast to display.

In addition, this class inherits all slots from its parent RowDotPlot, DotPlot, and Panel classes.

Examples

x <- VolcanoPlot()
x