mCSEATest: mCSEA core analysis

View source: R/mCSEATest.R

mCSEATestR Documentation

mCSEA core analysis

Description

Perform a methylated CpG sites enrichment analysis in predefined genomic regions

Usage

mCSEATest(rank, methData, pheno = NULL, column = 1,
  regionsTypes = c("promoters", "genes", "CGI"), customAnnotation = NULL,
  minCpGs = 5, nproc = 1, nperm = NULL, platform = "450k")

Arguments

rank

A named numeric vector with the ranking statistic of each CpG site

methData

A data frame or a matrix containing Illumina's CpG probes in rows and samples in columns. A SummarizedExperiment object can be used too

pheno

A data frame or a matrix containing samples in rows and covariates in columns. If NULL (default), pheno is extracted from the SummarizedExperiment object

column

The column name or number from pheno used to split the samples into groups (first column is used by default)

regionsTypes

A character or character vector indicating the predefined regions to be analyzed. NULL to skip this step and use customAnnotation.

customAnnotation

An optional list with the CpGs associated to each feature (default = NULL)

minCpGs

Minimum number of CpGs associated to a region. Regions below this threshold are not tested

nproc

Number of processors to use in GSEA step (default = 1)

nperm

(deprecated) Number of permutations to do in GSEA step in the previous implementation. Now, this parameter is ignored

platform

Platform used to get the methylation data (450k or EPIC)

Value

A list with the results of each of the analyzed regions. For each region type, a data frame with the results and a list with the probes associated to each region are generated. In addition, this list also contains the input methData, pheno and platform objects

Author(s)

Jordi Martorell Marugán, jordi.martorell@genyo.es

References

Subramanian, A. et al (2005). Gene set enrichment analysis: A knowledge-based approach for interpreting genome-wide expression profiles . PNAS 102, 15545-15550.

See Also

rankProbes, mCSEAPlot, mCSEAPlotGSEA

Examples

## Not run: 
library(mCSEAdata)
data(mcseadata)
myRank <- rankProbes(betaTest, phenoTest, refGroup = "Control")
set.seed(123)
myResults <- mCSEATest(myRank, betaTest, phenoTest,
regionsTypes = "promoters", platform = "EPIC")

## End(Not run)
data(precomputedmCSEA)
head(myResults[["promoters"]])
head(myResults[["promoters_association"]])

jordimartorell/mCSEA documentation built on Oct. 24, 2024, 4:44 a.m.