callBinary: callBinary

Description Usage Arguments Details Value Examples

View source: R/wrappers.R

Description

One of two main functions in the chromswitch package, this function detects a switch in chromatin state in one or more regions given ChIP-seq peak calls for one mark, executing the entire algorithm from preprocessing to evaluating the clustering results, using the binary strategy.

Usage

1
2
3
4
5
callBinary(query, metadata, peaks, filter = FALSE,
  filter_columns = NULL, filter_thresholds = NULL, reduce = TRUE,
  gap = 300, p = 0.4, n_features = FALSE, heatmap = FALSE,
  titles = NULL, outdir = NULL, optimal_clusters = TRUE,
  estimate_state = FALSE, test_condition = NULL, BPPARAM = bpparam())

Arguments

query

GRanges list containing one or more genomic regions of interest in which to call a switch. The output dataframe will contain one row per region in query.

metadata

A dataframe with at least two columns: "Sample" which stores the sample IDs, "Condition", which stores the biological condition labels of the samples

peaks

List of GRanges objects storing peak calls for each sample, where element names correspond to sample IDs

filter

(Optional) logical value, filter peaks based on thresholds on peak statistics? Default: FALSE. The filter step is described in filterPeaks.

filter_columns

If filter is TRUE, a chracter vector corresponding to names of columns in the peak metadata by which to filter peaks. If filter is FALSE, not used.

filter_thresholds

If filter is TRUE, a numeric vector corresponding to lower cutoffs applied to metadata columns in order to filter peaks. Provide one per column specified in filter_columns, in the same order. If filter is FALSE, not used.

reduce

(Optional) logical value, if TRUE, reduce gaps between nearby peaks in the same sample. See more at reducePeaks. Default: TRUE

gap

(Optional) If reduce is TRUE, numeric value, specifying the threshold distance for merging. Peaks in the same sample which are within this many bp of each other will be merged. Default: 300

p

Numeric value in [0, 1] giving the fraction of reciprocal overlap to require. Default: 0.4

n_features

(Optional) Logical value indicating whether to include a column "n_features" in the output storing the number of features in the feature matrix constructed for the region, which may be useful for understanding the behaviour of the binary strategy for constructing feature matrices. Default: FALSE

heatmap

(Optional) Logical value, plot the heatmap corresponding to the hierarchical clustering result? Default: FALSE

titles

(Optional) if heatmap is TRUE, a character vector of the same length as query, specifying the title to use when plotting each heatmap (e.g. a gene name), also reused as the prefix of the name of the file where the heatmap is saved. By default, the title is the genomic coordinates of the region in the form "chrN:start-end"

outdir

(Optional) if heatmap is TRUE, the name of the directory where heatmaps should be saved

optimal_clusters

(Optional) Logical value indicate whether to cluster samples into two groups, or to find the optimal clustering solution by choosing the set of clusters which maximizes the Average Silhouette width. Default: TRUE.

estimate_state

(Optional) Logical value indicating whether to include a column "state" in the output specifying the estimated chromatin state of a test condition. The state will be on of "ON", "OFF", or NA, where the latter results if a binary switch between the conditions is unclear. Default: FALSE.

test_condition

(Optional) If estimate_state is TRUE, string specifying one of the two biological condtions in metadata$Condition for which to estimate chromatin state.

BPPARAM

(Optional) instance of BiocParallel:BiocParallelParam used to determine the back-end used for parallel computations when performing the analysis on more than one region.

Details

This strategy constructs a sample-by-feature matrix to use as input for hierarchical clustering by first assembling the set of unique peaks observed in the region across samples. Then for each unique peak, we model the presence or absence of that peak in each sample, resulting in a binary feature matrix.

Value

Data frame with one row per region in query. Contains the coordinates of the region, the number of inferred clusters, the computed cluster validity statistics, and the cluster assignment for each sample.

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
samples <- c("E068", "E071", "E074", "E101", "E102", "E110")
bedfiles <- system.file("extdata", paste0(samples, ".H3K4me3.bed"),
package = "chromswitch")
Conditions <- c(rep("Brain", 3), rep("Other", 3))

metadata <- data.frame(Sample = samples,
    H3K4me3 = bedfiles,
    Condition = Conditions,
    stringsAsFactors = FALSE)

regions <- GRanges(seqnames = c("chr19", "chr19"),
    ranges = IRanges(start = c(54924104, 54874318),
                                end = c(54929104, 54877536)))

callBinary(query = regions, metadata = metadata, peaks = H3K4me3,
           BPPARAM = BiocParallel::SerialParam())

chromswitch documentation built on Nov. 8, 2020, 5:24 p.m.