computePvalues: Calculate P-values

computePvaluesR Documentation

Calculate P-values

Description

This function computes the P-values based on the fitted negative binomial model being an outlier for the given sample gene combination. It computes two matrices with the same dimension as the count matrix (samples x genes), which contain the corresponding P-values and adjusted P-values of every count. The adjusted P-values are computed across all genes per sample.

Usage

computePvalues(object, ...)

## S4 method for signature 'OutriderDataSet'
computePvalues(
  object,
  alternative = c("two.sided", "greater", "less"),
  method = "BY",
  subsets = NULL,
  BPPARAM = bpparam()
)

Arguments

object

An OutriderDataSet

...

additional params, currently not used.

alternative

Can be one of "two.sided", "greater" or "less" to specify the alternative hypothesis used to calculate the P-values, defaults to "two.sided"

method

Method used for multiple testing

subsets

A named list of named lists specifying any number of gene subsets (can differ per sample). For each subset, FDR correction will be limited to genes in the subset, and the FDR corrected pvalues stored as an assay in the ods object in addition to the transcriptome-wide FDR corrected pvalues. See the examples for how to use this argument.

BPPARAM

Can be used to parallelize the computation, defaults to bpparam()

Value

An OutriderDataSet object with computed nominal and adjusted P-values

See Also

p.adjust

Examples

ods <- makeExampleOutriderDataSet()

ods <- estimateSizeFactors(ods)
ods <- fit(ods)

ods <- computePvalues(ods)

assays(ods)[['pValue']][1:10,1:10]

# example of restricting FDR correction to subsets of genes of interest
genesOfInterest <- list("sample_1"=sample(rownames(ods), 3), 
                         "sample_2"=sample(rownames(ods), 8), 
                         "sample_6"=sample(rownames(ods), 5))
ods <- computePvalues(ods, subsets=list("exampleSubset"=genesOfInterest))
padj(ods, subsetName="exampleSubset")[1:10,1:10]
ods <- computePvalues(ods, 
                 subsets=list("anotherExampleSubset"=rownames(ods)[1:5]))
padj(ods, subsetName="anotherExampleSubset")[1:10,1:10]


gagneurlab/OUTRIDER documentation built on April 29, 2024, 2:22 a.m.