Find motif matches

matchMotifs(pwms, subject, ...)

# S4 method for PWMatrixList,DNAStringSet
matchMotifs(pwms, subject,
  genome = NULL, bg = c("subject", "genome", "even"), out = c("matches",
  "scores", "positions"), p.cutoff = 5e-05, w = 7, ranges = NULL)

# S4 method for PWMatrixList,character
matchMotifs(pwms, subject, genome = NULL,
  bg = c("subject", "genome", "even"), out = c("matches", "scores",
  "positions"), p.cutoff = 5e-05, w = 7, ranges = NULL)

# S4 method for PWMatrixList,DNAString
matchMotifs(pwms, subject, genome = NULL,
  bg = c("subject", "genome", "even"), out = c("matches", "scores",
  "positions"), p.cutoff = 5e-05, w = 7, ranges = NULL)

# S4 method for PWMatrixList,GenomicRanges
matchMotifs(pwms, subject,
  genome = GenomeInfoDb::genome(subject), bg = c("subject", "genome",
  "even"), out = c("matches", "scores", "positions"), p.cutoff = 5e-05,
  w = 7)

# S4 method for PWMatrixList,RangedSummarizedExperiment
matchMotifs(pwms, subject,
  genome = GenomeInfoDb::genome(subject), bg = c("subject", "genome",
  "even"), out = c("matches", "scores", "positions"), p.cutoff = 5e-05,
  w = 7)

# S4 method for PWMatrixList,BSgenomeViews
matchMotifs(pwms, subject,
  bg = c("subject", "genome", "even"), out = c("matches", "scores",
  "positions"), p.cutoff = 5e-05, w = 7)

# S4 method for PFMatrixList,ANY
matchMotifs(pwms, subject, ...)

# S4 method for PWMatrix,ANY
matchMotifs(pwms, subject, ...)

# S4 method for PFMatrix,ANY
matchMotifs(pwms, subject, ...)

Arguments

pwms
either PFMatrix, PFMatrixList, PWMatrix, PWMatrixList
subject
either GenomicRanges, DNAStringSet, DNAString, or character vector
...
additional arguments depending on inputs
genome
BSgenome object, DNAStringSet, or FaFile, or short string signifying genome build recognized by getBSgenome. Only required if subject is GenomicRanges or RangedSummarizedExperiment or if bg is set to "genome"
bg
background nucleotide frequencies. Default is to compute based on subject, i.e. the specific set of sequences being evaluated. See Details.
out
what to return? see return section
p.cutoff
p-value cutoff for returning motifs
w
parameter controlling size of window for filtration; default is 7
ranges
if subject is not GenomicRanges or RangedSummarizedExperiment, these ranges can be used to specify what ranges the input sequences correspond to. These ranges will be incorporated into the SummarizedExperiment output if out is "matches" or "scores" or will be used to give absolute positions of motifs if out is "positions"

Value

Either returns a SummarizedExperiment with a sparse matrix with values set to TRUE for a match (if out == 'matches'), a SummarizedExperiment with a matches matrix as well as matrices with the maximum motif score and total motif counts (if out == 'scores'), or a GenomicRangesList or a list of IRangesList with all the positions of matches (if out == 'positions')

Details

Background nucleotide frequencies can be set to "subject" to use the subject sequences or ranges for computing the nucleotide frequencies, "genome" for using the genomice frequencies (in which case a genome must be specified), "even" for using 0.25 for each base, or a numeric vector with A, C, G, and T frequencies.

Methods (by class)

  • pwms = PWMatrixList,subject = DNAStringSet: PWMatrixList/DNAStringSet

  • pwms = PWMatrixList,subject = character: PWMatrixList/character

  • pwms = PWMatrixList,subject = DNAString: PWMatrixList/DNAString

  • pwms = PWMatrixList,subject = GenomicRanges: PWMatrixList/GenomicRanges

  • pwms = PWMatrixList,subject = RangedSummarizedExperiment: PWMatrixList/RangedSummarizedExperiment

  • pwms = PWMatrixList,subject = BSgenomeViews: PWMatrixList/BSGenomeViews

  • pwms = PFMatrixList,subject = ANY: PFMatrixList/ANY

  • pwms = PWMatrix,subject = ANY: PWMatrix/ANY

  • pwms = PFMatrix,subject = ANY: PFMatrix/ANY

Examples

data(example_motifs, package = "motifmatchr") # Make a set of peaks peaks <- GenomicRanges::GRanges(seqnames = c("chr1","chr2","chr2"), ranges = IRanges::IRanges(start = c(76585873,42772928, 100183786), width = 500)) # Get motif matches for example motifs motif_ix <- matchMotifs(example_motifs, peaks, genome = "BSgenome.Hsapiens.UCSC.hg19")