motifScan: motifScan

motifScanR Documentation

motifScan

Description

function to find motif matches

Usage

motifScan(pwms, subject, ...)

## S4 method for signature 'PWMatrixList,DNAStringSet'
motifScan(
  pwms,
  subject,
  genome = NULL,
  bg = c("genome", "subject", "even"),
  out = c("matches", "scores", "positions"),
  p.cutoff = 1e-04,
  ranges = NULL,
  thread = 1,
  random.seed = NULL,
  cutoff.matrix.loc = "./",
  cutoff.matrix.name = NULL
)

## S4 method for signature 'PWMatrixList,character'
motifScan(
  pwms,
  subject,
  genome = NULL,
  bg = c("genome", "subject", "even"),
  out = c("matches", "scores", "positions"),
  p.cutoff = 1e-04,
  ranges = NULL,
  thread = 1,
  random.seed = NULL,
  cutoff.matrix.loc = "./",
  cutoff.matrix.name = NULL
)

## S4 method for signature 'PWMatrixList,DNAString'
motifScan(
  pwms,
  subject,
  genome = NULL,
  bg = c("genome", "subject", "even"),
  out = c("matches", "scores", "positions"),
  p.cutoff = 1e-04,
  ranges = NULL,
  thread = 1,
  random.seed = NULL,
  cutoff.matrix.loc = "./",
  cutoff.matrix.name = NULL
)

## S4 method for signature 'PWMatrixList,GenomicRanges'
motifScan(
  pwms,
  subject,
  genome = GenomeInfoDb::genome(subject),
  bg = c("genome", "subject", "even"),
  out = c("matches", "scores", "positions"),
  p.cutoff = 1e-04,
  thread = 1,
  random.seed = NULL,
  cutoff.matrix.loc = "./",
  cutoff.matrix.name = NULL
)

## S4 method for signature 'PWMatrixList,RangedSummarizedExperiment'
motifScan(
  pwms,
  subject,
  genome = GenomeInfoDb::genome(subject),
  bg = c("genome", "subject", "even"),
  out = c("matches", "scores", "positions"),
  p.cutoff = 1e-04,
  thread = 1,
  random.seed = NULL,
  cutoff.matrix.loc = "./",
  cutoff.matrix.name = NULL
)

## S4 method for signature 'PWMatrixList,BSgenomeViews'
motifScan(
  pwms,
  subject,
  bg = c("genome", "subject", "even"),
  out = c("matches", "scores", "positions"),
  p.cutoff = 1e-04,
  thread = 1,
  random.seed = NULL,
  cutoff.matrix.loc = "./",
  cutoff.matrix.name = NULL
)

## S4 method for signature 'PFMatrixList,ANY'
motifScan(pwms, subject, ...)

## S4 method for signature 'PWMatrix,ANY'
motifScan(pwms, subject, ...)

## S4 method for signature 'PFMatrix,ANY'
motifScan(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, 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". if the bg is set to "genome" and genome is FaFile, this function could only work on linux

bg

background nucleotide frequencies. Default is to compute based on genome, i.e. the specific genome being evaluated. See Details.

out

what to return? see value section

p.cutoff

p-value cutoff for returning motifs, should be one of 0.01, 0.001, 1e-04, 1e-05, 1e-06.(default=1e-04)

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"

thread

thread for running motifscan

random.seed

seed number for random program

cutoff.matrix.loc

the location of local motif score cutoff file, if the file is not present, motifscanR will generate by itself and save in the current working directory as './species_collect_cutoff_motifs_matrix.Rdata'(default), and the user could specify a specific file directory with 'save_path', or specify the cutoff matrix file by user himself with parameter cutoff.matrix.name

cutoff.matrix.name

the name of local motif score cutoff file, if the file is not present, motifscanR will generate by itself and save in the current working directory as './species_collect_cutoff_motifs_matrix.Rdata'(default), and the user could specify a specific file directory with 'save_path', then save as './[cutoff.matrix.name]_cutoff_motifs_matrix.Rdata, so user could use the next time if the pwms and the genome are same.

Details

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

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')

Methods (by class)

  • motifScan(pwms = PWMatrixList, subject = DNAStringSet): PWMatrixList/DNAStringSet

  • motifScan(pwms = PWMatrixList, subject = character): PWMatrixList/character

  • motifScan(pwms = PWMatrixList, subject = DNAString): PWMatrixList/DNAString

  • motifScan(pwms = PWMatrixList, subject = GenomicRanges): PWMatrixList/GenomicRanges

  • motifScan(pwms = PWMatrixList, subject = RangedSummarizedExperiment): PWMatrixList/RangedSummarizedExperiment

  • motifScan(pwms = PWMatrixList, subject = BSgenomeViews): PWMatrixList/BSGenomeViews

  • motifScan(pwms = PFMatrixList, subject = ANY): PFMatrixList/ANY

  • motifScan(pwms = PWMatrix, subject = ANY): PWMatrix/ANY

  • motifScan(pwms = PFMatrix, subject = ANY): PFMatrix/ANY

Examples

example_motifs <- getJasparMotifs(species = "Homo sapiens",
                                  collection = "CORE")

# Make a set of peaks
peaks <- GenomicRanges::GRanges(seqnames = c("chr1","chr2","chr2"),
                ranges = IRanges::IRanges(start = c(76585873,42772928,
                                          100183786),
                                          width = 500))

# Scan motif for example motifs
motif_ix <- motifScan(example_motifs, peaks, genome = "BSgenome.Hsapiens.UCSC.hg19")


LouisKwok-PICB/motifscanr documentation built on July 22, 2024, 6:36 a.m.