motifEnrichment: Enrichment of motif hits

Description Usage Arguments Details Value See Also Examples

View source: R/enrichmentTest.R

Description

This function determines whether a given motif is enriched in a given DNA sequences.

Usage

1
motifEnrichment(seqs, pfm, bg, singlestranded = FALSE, method = "compound")

Arguments

seqs

A DNAStringSet or DNAString object

pfm

An R matrix that represents a position frequency matrix

bg

A Background object

singlestranded

Boolean that indicates whether a single strand or both strands shall be scanned for motif hits. Default: singlestranded = FALSE.

method

String that defines whether to use the 'compound' Poisson approximation' or the 'combinatorial' model. Default: method='compound'.

Details

Enrichment is tested by comparing the observed number of motif hits against a theoretical distribution of the number of motif hits in random DNA sequences. Optionally, the theoretical distribution of the number of motif hits can be evaluated by either a 'compound Poisson model' or the 'combinatorial model'. Additionally, the enrichment test can be conducted with respect to scanning only the forward strand or both strands of the DNA sequences. The latter option is only available for the 'compound Poisson model'

Value

List that contains

pvalue

P-value for the enrichment test

fold

Fold-enrichment with respect to the expected number of hits

See Also

compoundPoissonDist, combinatorialDist

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
# Load sequences
seqfile = system.file("extdata", "seq.fasta", package = "motifcounter")
seqs = Biostrings::readDNAStringSet(seqfile)

# Load background
bg = readBackground(seqs, 1)

# Load motif
motiffile = system.file("extdata", "x31.tab", package = "motifcounter")
motif = t(as.matrix(read.table(motiffile)))

# 1 ) Motif enrichment test w.r.t. scanning a *single* DNA strand
# based on the 'Compound Poisson model'

result = motifEnrichment(seqs, motif, bg,
            singlestranded = TRUE, method = "compound")

# 2 ) Motif enrichment test w.r.t. scanning *both* DNA strand
# based on the 'Compound Poisson model'

result = motifEnrichment(seqs, motif, bg, method = "compound")

# 3 ) Motif enrichment test w.r.t. scanning *both* DNA strand
# based on the *combinatorial model*

result = motifEnrichment(seqs, motif, bg, singlestranded = FALSE,
            method = "combinatorial")

motifcounter documentation built on Nov. 8, 2020, 5:44 p.m.