View source: R/prepareBindingSites.R
prepareBindingSites | R Documentation |
Prepare binding sites by given position weight matrix and genome.
prepareBindingSites(
pwms,
genome,
seqlev = seqlevels(genome),
p.cutoff = 1e-05,
w = 7,
grange,
maximalBindingWidth = 40L,
mergeBindingSitesByPercentage = 0.8,
ignore.strand = TRUE
)
pwms |
either |
genome |
|
seqlev |
A character vector. Sequence levels to be searched. |
p.cutoff |
p-value cutoff for returning motifs; default is 1e-05 |
w |
parameter controlling size of window for filtration; default is 7 |
grange |
GRanges for motif search. If it is set, function will only search the binding site within the grange. Usually a peak list should be supplied. |
maximalBindingWidth |
A numeric vector(length=1). Maximal binding site width. Default is 40. |
mergeBindingSitesByPercentage |
A numeric vector (length=1). The percentage of overlapping region of binding sites to merge as one binding site. |
ignore.strand |
When set to TRUE, the strand information is ignored in the calculations. |
A GenomicRanges
with
all the positions of matches.
Jianhong Ou
library(TFBSTools)
motifs <- readRDS(system.file("extdata", "PWMatrixList.rds",
package="ATACseqTFEA"))
library(BSgenome.Drerio.UCSC.danRer10)
seqlev <- "chr1" #paste0("chr", 1:25)
mts <- prepareBindingSites(motifs, Drerio, seqlev,
grange=GRanges("chr1",
IRanges(5000, 100000)))
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.