biosCamera: An adapted and enhanced version of limma::camera

biosCameraR Documentation

An adapted and enhanced version of limma::camera

Description

An adapted and enhanced version of limma::camera

Usage

biosCamera(
  y,
  index,
  design = NULL,
  contrast = ncol(design),
  weights = NULL,
  geneLabels = NULL,
  use.ranks = FALSE,
  allow.neg.cor = FALSE,
  trend.var = FALSE,
  sort = FALSE,
  .fixed.inter.gene.cor = NULL,
  .approx.zscoreT = FALSE
)

Arguments

y

a numeric matrix of log-expression values or log-ratios of expression values, or any data object containing such a matrix. Rows correspond to probes and columns to samples. Any type of object that can be processed by getEAWP is acceptable.

index

an index vector or a list of index vectors. Can be any vector such that y[index,] of statistic[index] selects the rows corresponding to the test set. The list can be made using ids2indices.

design

Design matrix

contrast

contrast of the linear model coefficients for which the test is required. Can be an integer specifying a column of design, or else a numeric vector of same length as the number of columns of design.

weights

numeric matrix of observation weights of same size as y, or a numeric vector of array weights with length equal to ncol(y), or a numeric vector of gene weights with length equal to nrow(y).

geneLabels

Labels of the features in the input matrix.

use.ranks

do a rank-based test (TRUE) or a parametric test (FALSE)?

allow.neg.cor

should reduced variance inflation factors be allowed for negative correlations?

trend.var

logical, should an empirical Bayes trend be estimated? See eBayes for details.

sort

logical, should the results be sorted by p-value?

.fixed.inter.gene.cor

Numeric value, vector, or NULL/NA, advanced parameter corresponding to inter.gene.cor in the original implementation in limma. If set, gene-sets are set to have the fixed inter-gene correlation; the vector will be recycled to meet the correct length. If set as NULL/NA, correlations are estimated from each gene-set.

.approx.zscoreT

logical, advanced parameter only used for debugging purposes. If TRUE, the code is expected to return the exact same results as edgeR::camera (version 3.20.9), and maybe faster in execution.

The function was adapted from camera, with following improvments

  1. The output data.frame is more user-friendly

  2. The column 'FDR' is always present, even when only one gene-set was tested

  3. Scores are calculated, defined as log10(pValue)*I(directionality), where I(directionality) equals 1 if the directionality is Up and -1 if the directionality is Down

  4. Contributing genes and statistics are printed

Value

A data.frame with one row per set and the following columns:

GeneSet

Gene set name

NGenes

Number of genes in the set

Correlation

Estimated correlation

EffectSize

Estimated difference between the mean values of genes in the geneset and the background genes

Direction

Direction of set-wise regulation, Up or Down

Score

Gene-set enrichment score, defined as log10(pValue)*I(directionality), where I(directionality) equals 1 if the directionality is Up and -1 if the directionality is Down

ContribuingGenes

A character string, containing all genes labels of genes that are in the set and regulated in the same direction as the set-wise direction, and the respective statistic

Note

Since limma 3.29.6, the default setting of allow.neg.cor changes from TRUE to FALSE, and a new parameter, inter.gene.cor, is added with the default value of 0.01, namely a prior inter-gene correlation is set for all gene sets. Currently, biosCamera does not have the parameter inter.gene.cor, but allow.neg.cor is set by default to FALSE to be consistent with the latest camera function.

Examples

y <- matrix(rnorm(1000*6),1000,6)
design <- cbind(Intercept=1,Group=c(0,0,0,1,1,1))
# First set of 20 genes are genuinely deferentially expressed 
index1 <- 1:20
y[index1,4:6] <- y[index1,4:6]+1
# The second set of 20 genes are not
index2 <- 21:40
biosCamera(y, index1, design) 
biosCamera(y, index2, design)
biosCamera(y, list(index1, index2), design)

# compare with the output of camera: columns 'GeneSet', 'Score',
# 'ContributingGenes' are missing, and in case \code{inter.gene.cor} is (as
# default) set to a numeric value, the column 'Correlation' is also missing

limmaDefOut <- limma::camera(y, index1, design)
limmaCorDefOut <-
    limma::camera(y, index1, design, inter.gene.cor=NA)

## Not run:  
  # when \code{.approx.zscoreT=TRUE},  PValue reported by
  # \code{limma::camera(inter.gene.cor=NA)} and \code{ribiosGSEA::biosCamera}
  # should equal 
  biosCorOut <- biosCamera(y, index1, design, .approx.zscoreT=TRUE)

  # when \code{.fixed.inter.gene.cor=0.01} and \code{.approx.zscoreT=TRUE},
  # PValue reported by \code{limma::camera} and \code{ribiosGSEA::biosCamera}
  # should equal 
  biosFixCorOut <- biosCamera(y, index1, design,
      .fixed.inter.gene.cor=0.01, .approx.zscoreT=TRUE)
  testthat::expect_equal(biosFixCorOut$PValue, limmaDefOut$PValue)
  testthat::expect_equal(biosCorOut$PValue, limmaCorDefOut$PValue)

## End(Not run)


bedapub/ribiosGSEA documentation built on March 30, 2023, 3:26 p.m.