get_tf_in_region: Get human TFs for regions by either scanning it with...

Description Usage Arguments Value Examples

View source: R/get_tf_in_region.R

Description

Given a genomic region, this function maps TF in regions using two methods: 1) using motifmatchr nd JASPAR2020 to scan the region for 554 human transcription factors binding sites. There is also an option (argument window.size) to extend the scanning region before performing the search, which by default is 0 (do not extend). 2) Using user input TF chip-seq to check for overlaps between region and TF peaks.

Usage

1
2
3
4
5
6
7
8
9
get_tf_in_region(
  region,
  window.size = 0,
  genome = c("hg19", "hg38"),
  p.cutoff = 1e-08,
  cores = 1,
  TF.peaks.gr = NULL,
  verbose = FALSE
)

Arguments

region

A vector of region names or GRanges object with the DNA methylation regions to be scanned for the motifs

window.size

Integer value to extend the regions. For example, a value of 50 will extend 25 bp upstream and 25 bp downstream the region. The default is not to increase the scanned region.

genome

Human genome of reference "hg38" or "hg19".

p.cutoff

motifmatchr p.cutoff. Default 1e-8.

cores

Number of CPU cores to be used. Default 1.

TF.peaks.gr

A granges with TF peaks to be overlaped with input region Metadata column expected "id" with TF name. Default NULL. Note that Remap catalog can be used as shown in the examples.

verbose

A logical argument indicating if messages output should be provided.

Value

A data frame with the following information: regionID, TF symbol, TF ensembl ID

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
 regions.names <- c("chr3:189631389-189632889","chr4:43162098-43163498")
 region.tf <- get_tf_in_region(
                 region = regions.names,
                 genome = "hg38"
 )

## Not run: 
   library(ReMapEnrich)
   demo.dir <- "~/ReMapEnrich_demo"
   dir.create(demo.dir, showWarnings = FALSE, recursive = TRUE)
   # Use the function DowloadRemapCatalog
   remapCatalog2018hg38 <- downloadRemapCatalog(demo.dir, assembly = "hg38")
   # Load the ReMap catalogue and convert it to Genomic Ranges
   remapCatalog <- bedToGranges(remapCatalog2018hg38)
   regions.names <- c("chr3:189631389-189632889","chr4:43162098-43163498")
   region.tf.remap <- get_tf_in_region(
                   region = regions.names,
                   genome = "hg38",
                   TF.peaks.gr = remapCatalog
   )
 
## End(Not run)

MethReg documentation built on Nov. 8, 2020, 8:01 p.m.