deseq_polar: Convert DESeq2 objects to a volcano3d object

View source: R/deseq_polar.R

deseq_polarR Documentation

Convert DESeq2 objects to a volcano3d object

Description

This function is used instead of polar_coords if you have raw RNA-Seq count data. It takes 2 DESeqDataSet objects, extracts statistical results and converts the results to a 'volc3d' object, which can be directly plotted.

Usage

deseq_polar(
  object,
  objectLRT,
  contrast = NULL,
  data = NULL,
  pcutoff = 0.05,
  padj.method = "BH",
  filter_pairwise = TRUE,
  ...
)

Arguments

object

An object of class 'DESeqDataSet' with the full design formula. The function DESeq needs to have been run.

objectLRT

An object of class 'DESeqDataSet' with the reduced design formula. The function DESeq needs to have been run on this object with argument test="LRT".

contrast

Character value specifying column within the metadata stored in the DESeq2 dataset objects is the outcome variable. This column must contain a factor with 3 levels. If not set, the function will select the last term in the design formula of object as per DESeq2 convention.

data

Optional matrix containing gene expression data. If not supplied, the function will pull the expression data from within the DESeq2 object using the DESeq2 function assay(). NOTE: for consistency with gene expression datasets, genes are in rows.

pcutoff

Cut-off for p-value significance

padj.method

Can be any method available in p.adjust or "qvalue". The option "none" is a pass-through.

filter_pairwise

Logical whether adjusted p-value pairwise statistical tests are only conducted on genes which reach significant adjusted p-value cut-off on the group likelihood ratio test

...

Optional arguments passed to polar_coords

Value

Calls polar_coords to return an S4 'volc3d' object

See Also

polar_coords, voom_polar, DESeq in the DESeq2 package

Examples



  library(DESeq2)

  counts <- matrix(rnbinom(n=1500, mu=100, size=1/0.5), ncol=15)
  cond <- factor(rep(1:3, each=5), labels = c('A', 'B', 'C'))

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

  # standard analysis
  dds <- DESeq(dds)

  # Likelihood ratio test
  ddsLRT <- DESeq(dds, test="LRT", reduced= ~ 1)

  polar <- deseq_polar(dds, ddsLRT, "cond")
  volcano3D(polar)
  radial_ggplot(polar)



volcano3D documentation built on May 31, 2023, 5:31 p.m.