getDifferentialAccessibleTiles {MOCHA}R Documentation

getDifferentialAccessibleTiles

Description

getDifferentialAccessibleTiles allows you to determine whether regions of chromatin are differentially accessible between groups by conducting a test

Usage

getDifferentialAccessibleTiles(
  SampleTileObj,
  cellPopulation,
  groupColumn,
  foreground,
  background,
  signalThreshold = 12,
  minZeroDiff = 0.5,
  fdrToDisplay = 0.2,
  outputGRanges = TRUE,
  numCores = 1,
  verbose = FALSE
)

Arguments

SampleTileObj

The SummarizedExperiment object output from getSampleTileMatrix

cellPopulation

A string denoting the cell population of interest

groupColumn

The column containing sample group labels

foreground

The foreground group of samples for differential comparison

background

The background group of samples for differential comparison

signalThreshold

Minimum median intensity required to keep tiles for differential testing to increase statistical power in small sample cohorts. Default is 12.

minZeroDiff

Minimum difference in average dropout rates across groups require to keep tiles for differential testing. Default is 0.5 (50%).

fdrToDisplay

False-discovery rate used only for standard output messaging. Default is 0.2.

outputGRanges

Outputs a GRanges if TRUE and a data.frame if FALSE. Default is TRUE.

numCores

The number of cores to use with multiprocessing. Default is 1.

verbose

Set TRUE to display additional messages. Default is FALSE.

Value

full_results The differential accessibility results as a GRanges or matrix data.frame depending on the flag 'outputGRanges'.

Examples

## Not run: 
cellPopulation <- "MAIT"
foreground <- "Positive"
background <- "Negative"
# Standard output will display the number of tiles found below a false-discovery rate threshold.
# This parameter does not filter results and only affects the aforementioned message.
fdrToDisplay <- 0.2
# Choose to output a GRanges or data.frame.
# Default is TRUE
outputGRanges <- TRUE
# SampleTileMatrices is the output of MOCHA::getSampleTileMatrix
differentials <- MOCHA::getDifferentialAccessibleTiles(
  SampleTileObj = SampleTileMatrices,
  cellPopulation = cellPopulation,
  groupColumn = groupColumn,
  foreground = foreground,
  background = background,
  fdrToDisplay = fdrToDisplay,
  outputGRanges = outputGRanges,
  numCores = numCores
)

## End(Not run)

[Package MOCHA version 1.1.0 Index]