Title: | Identification of enriched motif pairs from chromatin interaction data |
---|---|
Description: | Identifies motifs that are significantly co-enriched from enhancer-promoter interaction data. While enhancer-promoter annotation is commonly used to define groups of interaction anchors, spatzie also supports co-enrichment analysis between preprocessed interaction anchors. Supports BEDPE interaction data derived from genome-wide assays such as HiC, ChIA-PET, and HiChIP. Can also be used to look for differentially enriched motif pairs between two interaction experiments. |
Authors: | Jennifer Hammelman [aut, cre, cph] , Konstantin Krismer [aut] , David Gifford [ths, cph] |
Maintainer: | Jennifer Hammelman <[email protected]> |
License: | GPL-3 |
Version: | 1.13.0 |
Built: | 2024-12-14 04:10:37 UTC |
Source: | https://github.com/bioc/spatzie |
Determine whether motifs between paired bed regions have a statistically significant relationship. Options for significance are motif score correlation, motif count correlation, or hypergeometric motif co-occurrence.
anchor_pair_enrich(interaction_data, method = c("count", "score", "match"))
anchor_pair_enrich(interaction_data, method = c("count", "score", "match"))
interaction_data |
an interactionData object of paired genomic regions |
||||||
method |
method for co-occurrence, valid options include:
|
an interactionData object where obj$pair_motif_enrich
contains
the p-values for significance of seeing a higher co-occurrence than
what we get by chance.
We assume motif scores follow a normal distribution and are independent between enhancers and promoters. We can therefore compute how correlated scores of any two transcription factor motifs are between enhancer and promoter regions using Pearson's product-moment correlation coefficient:
,
where the input vectors and
from
above are transformed to vectors
and
by replacing the set of scores with the
maximum score for each region:
is then the maximum motif score of motif
in the
promoter region of interaction
,
is the maximum
motif score of motif
in the enhancer region of interaction
,
and
and
are the sample means.
Significance is then computed by transforming the correlation coefficient
to test statistic
, which is Student
-distributed
with
degrees of freedom.
All p-values are calculated as one-tailed p-values of the probability
that scores are greater than or equal to .
Instead of calculating the correlation of motif scores directly, the
count-based correlation metric first tallies the number of instances of a
given motif within an enhancer or a promoter region, which are defined as
all positions in those regions with motif score p-values of less than
. Formally, the input vectors
and
are transformed to vectors
and
by replacing the set of scores with the cardinality of the set:
And analogous for . Finally, the correlation
coefficient
between
and
and its associated significance are
calculated as described above.
Instance co-occurrence uses the presence or absence of a motif within an
enhancer or promoter to determine a statistically significant association,
thus and
are defined by:
Instance co-occurrence is computed using the hypergeometric test:
where is the number of interactions that contain a match for
motif
in the promoter and motif
in the enhancer,
is the number of promoters that contain motif
(
),
is the number of
enhancers that contain motif
(
), and
is the total
number of interactions, which is equal to the number of promoters and to
the number of enhancers.
Jennifer Hammelman
Konstantin Krismer
## Not run: genome_id <- "BSgenome.Mmusculus.UCSC.mm9" if (!(genome_id %in% rownames(utils::installed.packages()))) { BiocManager::install(genome_id, update = FALSE, ask = FALSE) } genome <- BSgenome::getBSgenome(genome_id) motifs_file <- system.file("extdata/motifs_subset.txt.gz", package = "spatzie") motifs <- TFBSTools::readJASPARMatrix(motifs_file, matrixClass = "PFM") yy1_pd_interaction <- scan_motifs(spatzie::interactions_yy1, motifs, genome) yy1_pd_interaction <- filter_motifs(yy1_pd_interaction, 0.4) yy1_pd_count_corr <- anchor_pair_enrich(yy1_pd_interaction, method = "count") ## End(Not run) res <- anchor_pair_enrich(spatzie::scan_interactions_example_filtered, method = "score")
## Not run: genome_id <- "BSgenome.Mmusculus.UCSC.mm9" if (!(genome_id %in% rownames(utils::installed.packages()))) { BiocManager::install(genome_id, update = FALSE, ask = FALSE) } genome <- BSgenome::getBSgenome(genome_id) motifs_file <- system.file("extdata/motifs_subset.txt.gz", package = "spatzie") motifs <- TFBSTools::readJASPARMatrix(motifs_file, matrixClass = "PFM") yy1_pd_interaction <- scan_motifs(spatzie::interactions_yy1, motifs, genome) yy1_pd_interaction <- filter_motifs(yy1_pd_interaction, 0.4) yy1_pd_count_corr <- anchor_pair_enrich(yy1_pd_interaction, method = "count") ## End(Not run) res <- anchor_pair_enrich(spatzie::scan_interactions_example_filtered, method = "score")
This object contains genomic interactions obtained by mouse YY1 ChIA-PET scanned for mouse transcription factor motifs, filtered for motifs present in at least 10 interactions with count correlation. It serves as unit test data.
data(anchor_pair_example_count)
data(anchor_pair_example_count)
An interactionData object
This object contains genomic interactions obtained by mouse YY1 ChIA-PET scanned for mouse transcription factor motifs, filtered for motifs present in at least 10 interactions with using the hypergeometric test. It serves as unit test data.
data(anchor_pair_example_match)
data(anchor_pair_example_match)
A interactionData object
This object contains genomic interactions obtained by mouse YY1 ChIA-PET scanned for mouse transcription factor motifs, filtered for motifs present in at least 10 interactions with score correlation. It serves as unit test data.
data(anchor_pair_example_score)
data(anchor_pair_example_score)
An interactionData object
Compute the log-likelihood ratio that a motif pair is differential between two interaction datasets. Note that motif pair significance should have been computed using the same method for both datasets.
compare_motif_pairs( interaction_data1, interaction_data2, differential_p = 0.05 )
compare_motif_pairs( interaction_data1, interaction_data2, differential_p = 0.05 )
interaction_data1 |
an interactionData object of paired genomic regions that has been scanned for significant motif:motif interactions |
interaction_data2 |
an interactionData object of paired genomic regions that has been scanned for significant motif:motif interactions |
differential_p |
threshold for significance of differential p-value |
a matrix of the log likelihood ratio of motif pairs that are significantly differential between two interactionData sets
Jennifer Hammelman
pheatmap::pheatmap(compare_motif_pairs(spatzie::int_data_k562, spatzie::int_data_mslcl, 5e-06), fontsize = 6)
pheatmap::pheatmap(compare_motif_pairs(spatzie::int_data_k562, spatzie::int_data_mslcl, 5e-06), fontsize = 6)
This is a matrix containing example result from compare_motif_pairs. It serves as unit test data.
data(compare_pairs_example)
data(compare_pairs_example)
A matrix
Select a subset of motifs that are in at least a threshold fraction of regions. Motif subsets are selected separately for anchor one and anchor two regions.
filter_motifs(interaction_data, threshold)
filter_motifs(interaction_data, threshold)
interaction_data |
an interactionData object of paired genomic regions |
threshold |
fraction of interactions that should contain a motif for a motif to be considered |
an interactionData object where obj$anchor1_motif_indices
and obj$anchor2_motif_indices
have been filtered to motifs that are
present in a threshold fraction of interactions
Jennifer Hammelman
## Not run: genome_id <- "BSgenome.Mmusculus.UCSC.mm9" if (!(genome_id %in% rownames(utils::installed.packages()))) { BiocManager::install(genome_id, update = FALSE, ask = FALSE) } genome <- BSgenome::getBSgenome(genome_id) motifs_file <- system.file("extdata/motifs_subset.txt.gz", package = "spatzie") motifs <- TFBSTools::readJASPARMatrix(motifs_file, matrixClass = "PFM") yy1_pd_interaction <- scan_motifs(spatzie::interactions_yy1, motifs, genome) yy1_pd_interaction <- filter_motifs(yy1_pd_interaction, 0.4) ## End(Not run) res <- filter_motifs(spatzie::scan_interactions_example, threshold = 0.1)
## Not run: genome_id <- "BSgenome.Mmusculus.UCSC.mm9" if (!(genome_id %in% rownames(utils::installed.packages()))) { BiocManager::install(genome_id, update = FALSE, ask = FALSE) } genome <- BSgenome::getBSgenome(genome_id) motifs_file <- system.file("extdata/motifs_subset.txt.gz", package = "spatzie") motifs <- TFBSTools::readJASPARMatrix(motifs_file, matrixClass = "PFM") yy1_pd_interaction <- scan_motifs(spatzie::interactions_yy1, motifs, genome) yy1_pd_interaction <- filter_motifs(yy1_pd_interaction, 0.4) ## End(Not run) res <- filter_motifs(spatzie::scan_interactions_example, threshold = 0.1)
Multiple hypothesis correction applied to filter for significant motif interactions.
filter_pair_motifs(interaction_data, method = "fdr", threshold = 0.05)
filter_pair_motifs(interaction_data, method = "fdr", threshold = 0.05)
interaction_data |
an interactionData object of paired genomic regions |
method |
statistical method for multiple hypothesis correction,
defaults to Benjamini-Hochberg ( |
threshold |
p-value threshold for significance cut-off |
an interactionData object where obj$pair_motif_enrich
contains
multiple hypothesis corrected p-values for significance of seeing a
higher co-occurrence than what we get by chance and
obj$pair_motif_enrich_sig
contains only motifs that have at least one
significant interaction.
Jennifer Hammelman
## Not run: genome_id <- "BSgenome.Mmusculus.UCSC.mm9" if (!(genome_id %in% rownames(utils::installed.packages()))) { BiocManager::install(genome_id, update = FALSE, ask = FALSE) } genome <- BSgenome::getBSgenome(genome_id) motifs_file <- system.file("extdata/motifs_subset.txt.gz", package = "spatzie") motifs <- TFBSTools::readJASPARMatrix(motifs_file, matrixClass = "PFM") yy1_pd_interaction <- scan_motifs(spatzie::interactions_yy1, motifs, genome) yy1_pd_interaction <- filter_motifs(yy1_pd_interaction, 0.4) yy1_pd_score_corr <- anchor_pair_enrich(yy1_pd_interaction, method = "score") yy1_pd_score_corr_adj <- filter_pair_motifs(yy1_pd_score_corr) ## End(Not run) res <- filter_pair_motifs(spatzie::anchor_pair_example_count, threshold = 0.5)
## Not run: genome_id <- "BSgenome.Mmusculus.UCSC.mm9" if (!(genome_id %in% rownames(utils::installed.packages()))) { BiocManager::install(genome_id, update = FALSE, ask = FALSE) } genome <- BSgenome::getBSgenome(genome_id) motifs_file <- system.file("extdata/motifs_subset.txt.gz", package = "spatzie") motifs <- TFBSTools::readJASPARMatrix(motifs_file, matrixClass = "PFM") yy1_pd_interaction <- scan_motifs(spatzie::interactions_yy1, motifs, genome) yy1_pd_interaction <- filter_motifs(yy1_pd_interaction, 0.4) yy1_pd_score_corr <- anchor_pair_enrich(yy1_pd_interaction, method = "score") yy1_pd_score_corr_adj <- filter_pair_motifs(yy1_pd_score_corr) ## End(Not run) res <- filter_pair_motifs(spatzie::anchor_pair_example_count, threshold = 0.5)
This object contains genomic interactions obtained by mouse YY1 ChIA-PET scanned for mouse transcription factor motifs, filtered for motifs present in at least 10 interactions with score correlation, and filtered for pairs with p < 0.5. It serves as unit test data.
data(filter_pairs_example)
data(filter_pairs_example)
An interactionData object
Identifies co-enriched pairs of motifs in enhancer-promoter interactions selected from a data frame of general genomic interactions.
If identify_ep
: Promoters and enhancers are identified
using genomic annotations, where anchors close to promoter annotations
(within 2500 base pairs) are considered promoters and all other anchors are
considered gene-distal enhancers. Only interactions in
int_raw_data
between promoters and enhancers are used for motif
co-enrichment analysis.
If !identify_ep
: Instead of automatically identifying
promoters and enhancers based on genomic annotations, all interactions
in int_raw_data
must be preprocessed in a way that anchor 1 contains
promoters and anchor 2 contains enhancers. Motif
co-enrichment analysis is performed under this assumption.
Calls functions scan_motifs
, filter_motifs
,
and anchor_pair_enrich
internally.
find_ep_coenrichment( int_raw_data, motifs_file, motifs_file_matrix_format = c("pfm", "ppm", "pwm"), genome_id = c("hg38", "hg19", "mm9", "mm10"), identify_ep = TRUE, cooccurrence_method = c("count", "score", "match"), filter_threshold = 0.4 )
find_ep_coenrichment( int_raw_data, motifs_file, motifs_file_matrix_format = c("pfm", "ppm", "pwm"), genome_id = c("hg38", "hg19", "mm9", "mm10"), identify_ep = TRUE, cooccurrence_method = c("count", "score", "match"), filter_threshold = 0.4 )
int_raw_data |
a
|
||||||||||||
motifs_file |
JASPAR format matrix file containing multiple motifs to scan for, gz-zipped files allowed |
||||||||||||
motifs_file_matrix_format |
type of position-specific scoring matrices
in
|
||||||||||||
genome_id |
ID of genome assembly interactions in |
||||||||||||
identify_ep |
logical, set |
||||||||||||
cooccurrence_method |
method for co-occurrence, valid options include:
See |
||||||||||||
filter_threshold |
fraction of interactions that should contain a
motif for a motif to be considered, see |
a list with the following items:
int_data |
GenomicInteractions object;
promoter-enhancer interactions |
int_data_motifs : |
interactionData object; return value of
scan_motifs
|
filtered_int_data_motifs : |
interactionData object;
return value of filter_motifs
|
annotation_pie_chart : |
ggplot2 plot; return value of
plotInteractionAnnotations
|
motif_cooccurrence : |
interactionData object; return
value of anchor_pair_enrich
|
Jennifer Hammelman
Konstantin Krismer
## Not run: interactions_file <- system.file("extdata/yy1_interactions.bedpe.gz", package = "spatzie") motifs_file <- system.file("extdata/motifs_subset.txt.gz", package = "spatzie") df <- read.table(gzfile(interactions_file), header = TRUE, sep = "\t") res <- find_ep_coenrichment(df, motifs_file, motifs_file_matrix_format = "pfm", genome_id = "mm10") ## End(Not run)
## Not run: interactions_file <- system.file("extdata/yy1_interactions.bedpe.gz", package = "spatzie") motifs_file <- system.file("extdata/motifs_subset.txt.gz", package = "spatzie") df <- read.table(gzfile(interactions_file), header = TRUE, sep = "\t") res <- find_ep_coenrichment(df, motifs_file, motifs_file_matrix_format = "pfm", genome_id = "mm10") ## End(Not run)
Select interactions that contain anchor1_motif within anchor 1 and anchor2_motif within anchor 2.
get_specific_interactions( interaction_data, anchor1_motif = NULL, anchor2_motif = NULL )
get_specific_interactions( interaction_data, anchor1_motif = NULL, anchor2_motif = NULL )
interaction_data |
an interactionData object of paired genomic regions |
anchor1_motif |
Motif name from |
anchor2_motif |
Motif name from |
a GenomicInteractions
object
containing a subset subset of interactions that contain an instance of
anchor1_motif
in anchor 1 and anchor2_motif
in anchor 2
Jennifer Hammelman
## Not run: genome_id <- "BSgenome.Mmusculus.UCSC.mm9" if (!(genome_id %in% rownames(utils::installed.packages()))) { BiocManager::install(genome_id, update = FALSE, ask = FALSE) } genome <- BSgenome::getBSgenome(genome_id) motifs_file <- system.file("extdata/motifs_subset.txt.gz", package = "spatzie") motifs <- TFBSTools::readJASPARMatrix(motifs_file, matrixClass = "PFM") yy1_pd_interaction <- scan_motifs(spatzie::interactions_yy1, motifs, genome) yy1_pd_interaction <- filter_motifs(yy1_pd_interaction, 0.4) yy1_pd_count_corr <- anchor_pair_enrich(yy1_pd_interaction, method = "score") yy1_yy1_interactions <- get_specific_interactions( yy1_pd_interaction, anchor1_motif = "YY1", anchor2_motif = "YY1") ## End(Not run) res <- get_specific_interactions(spatzie::int_data_yy1, anchor1_motif = "YY1", anchor2_motif = "YY1")
## Not run: genome_id <- "BSgenome.Mmusculus.UCSC.mm9" if (!(genome_id %in% rownames(utils::installed.packages()))) { BiocManager::install(genome_id, update = FALSE, ask = FALSE) } genome <- BSgenome::getBSgenome(genome_id) motifs_file <- system.file("extdata/motifs_subset.txt.gz", package = "spatzie") motifs <- TFBSTools::readJASPARMatrix(motifs_file, matrixClass = "PFM") yy1_pd_interaction <- scan_motifs(spatzie::interactions_yy1, motifs, genome) yy1_pd_interaction <- filter_motifs(yy1_pd_interaction, 0.4) yy1_pd_count_corr <- anchor_pair_enrich(yy1_pd_interaction, method = "score") yy1_yy1_interactions <- get_specific_interactions( yy1_pd_interaction, anchor1_motif = "YY1", anchor2_motif = "YY1") ## End(Not run) res <- get_specific_interactions(spatzie::int_data_yy1, anchor1_motif = "YY1", anchor2_motif = "YY1")
This object contains genomic interactions obtained by human RAD21 ChIA-PET from K562 cells and serves as unit test data.
data(int_data_k562)
data(int_data_k562)
An interactionData
object
This object contains genomic interactions obtained by human RAD21 ChIA-PET from MSLCL cells and serves as unit test data.
data(int_data_mslcl)
data(int_data_mslcl)
An interactionData
object
This object contains genomic interactions obtained by mouse YY1 ChIA-PET and serves as example and unit test data.
data(int_data_yy1)
data(int_data_yy1)
An interactionData
object
This object contains genomic interactions obtained by mouse YY1 ChIA-PET and serves as example and unit test data. The same data set is used in the vignette.
data(interactions_yy1)
data(interactions_yy1)
A GenomicInteractions
object
This is a GenomicInteractions
object
containing proccessed results from YY1
ChIA-PET of interactions that contain a YY1 motif in the enhancer
(anchor 2) region. It serves as unit test data.
data(interactions_yy1_enhancer)
data(interactions_yy1_enhancer)
A GenomicInteractions
object
This is a GenomicInteractions
object
containing proccessed results from YY1
ChIA-PET of interactions that contain a YY1 motif in the promoter
(anchor 1) region and a YY1 motif in the enhancer (anchor 2) region.
It serves as unit test data.
data(interactions_yy1_ep)
data(interactions_yy1_ep)
A GenomicInteractions
object
This is a GenomicInteractions
object
containing proccessed results from YY1
ChIA-PET of interactions that contain a YY1 motif in the promoter
(anchor 1) region. It serves as unit test data.
data(interactions_yy1_promoter)
data(interactions_yy1_promoter)
A GenomicInteractions
object
Plots a histogram of motif values (either counts, instances, or scores) for anchor 1 and anchor 2 regions.
plot_motif_occurrence( interaction_data, method = c("counts", "instances", "scores") )
plot_motif_occurrence( interaction_data, method = c("counts", "instances", "scores") )
interaction_data |
an interactionData object of paired genomic regions |
method |
way to interpret motif matching for each anchor region as "counts" number of motifs per region, "instances" motif present or absent each region, or "scores" maximum motif PWM match score for each region |
plot containing histogram for each anchor
Jennifer Hammelman
## Not run: genome_id <- "BSgenome.Mmusculus.UCSC.mm9" if (!(genome_id %in% rownames(utils::installed.packages()))) { BiocManager::install(genome_id, update = FALSE, ask = FALSE) } genome <- BSgenome::getBSgenome(genome_id) motifs_file <- system.file("extdata/motifs_subset.txt.gz", package = "spatzie") motifs <- TFBSTools::readJASPARMatrix(motifs_file, matrixClass = "PFM") yy1_pd_interaction <- scan_motifs(spatzie::interactions_yy1, motifs, genome) yy1_pd_interaction <- filter_motifs(yy1_pd_interaction, 0.4) plot_motif_occurrence(yy1_pd_interaction,"counts") ## End(Not run) plot_motif_occurrence(spatzie::anchor_pair_example_score)
## Not run: genome_id <- "BSgenome.Mmusculus.UCSC.mm9" if (!(genome_id %in% rownames(utils::installed.packages()))) { BiocManager::install(genome_id, update = FALSE, ask = FALSE) } genome <- BSgenome::getBSgenome(genome_id) motifs_file <- system.file("extdata/motifs_subset.txt.gz", package = "spatzie") motifs <- TFBSTools::readJASPARMatrix(motifs_file, matrixClass = "PFM") yy1_pd_interaction <- scan_motifs(spatzie::interactions_yy1, motifs, genome) yy1_pd_interaction <- filter_motifs(yy1_pd_interaction, 0.4) plot_motif_occurrence(yy1_pd_interaction,"counts") ## End(Not run) plot_motif_occurrence(spatzie::anchor_pair_example_score)
This object contains genomic interactions obtained by mouse YY1 ChIA-PET scanned for mouse transcription factor motifs and serves as unit test data.
data(scan_interactions_example)
data(scan_interactions_example)
An interactionData object
This object contains genomic interactions obtained by mouse YY1 ChIA-PET scanned for mouse transcription factor motifs and filtered for motifs present in at least 10
data(scan_interactions_example_filtered)
data(scan_interactions_example_filtered)
An interactionData object
Uses motifmatchR to scan interaction regions for given motifs.
scan_motifs(int_data, motifs, genome)
scan_motifs(int_data, motifs, genome)
int_data |
a |
motifs |
a TFBS tools matrix of DNA binding motifs |
genome |
BSgenome object or DNAStringSet object, must match chromosomes from interaction data file |
an interaction data object where obj$anchor1_motifs
and
obj$anchor2_motifs
contain information about the scores and matches
to motifs from anchor one and anchor two of interaction data genomic regions
Jennifer Hammelman
## Not run: genome_id <- "BSgenome.Mmusculus.UCSC.mm9" if (!(genome_id %in% rownames(utils::installed.packages()))) { BiocManager::install(genome_id, update = FALSE, ask = FALSE) } genome <- BSgenome::getBSgenome(genome_id) motifs_file <- system.file("extdata/motifs_subset.txt.gz", package = "spatzie") motifs <- TFBSTools::readJASPARMatrix(motifs_file, matrixClass = "PFM") yy1_pd_interaction <- scan_motifs(spatzie::interactions_yy1, motifs, genome) ## End(Not run) motifs_file <- system.file("extdata/motifs_subset.txt.gz", package = "spatzie") motifs <- TFBSTools::readJASPARMatrix(motifs_file, matrixClass = "PFM") left <- GenomicRanges::GRanges( seqnames = c("chr1", "chr1", "chr1"), ranges = IRanges::IRanges(start = c(1, 15, 20), end = c(10, 35, 31))) right <- GenomicRanges::GRanges( seqnames = c("chr1", "chr2", "chr2"), ranges = IRanges::IRanges(start = c(17, 47, 41), end = c(28, 54, 53))) test_interactions <- GenomicInteractions::GenomicInteractions(left, right) # toy DNAStringSet to replace BSgenome object seqs <- c("chr1" = "CCACTAGCCACGCGTCACTGGTTAGCGTGATTGAAACTAAATCGTATGAAAATCC", "chr2" = "CTACAAACTAGGAATTTAGGCAAACCTGTGTTAAAATCTTAGCTCATTCATTAAT") toy_genome <- Biostrings::DNAStringSet(seqs, use.names = TRUE) res <- scan_motifs(test_interactions, motifs, toy_genome)
## Not run: genome_id <- "BSgenome.Mmusculus.UCSC.mm9" if (!(genome_id %in% rownames(utils::installed.packages()))) { BiocManager::install(genome_id, update = FALSE, ask = FALSE) } genome <- BSgenome::getBSgenome(genome_id) motifs_file <- system.file("extdata/motifs_subset.txt.gz", package = "spatzie") motifs <- TFBSTools::readJASPARMatrix(motifs_file, matrixClass = "PFM") yy1_pd_interaction <- scan_motifs(spatzie::interactions_yy1, motifs, genome) ## End(Not run) motifs_file <- system.file("extdata/motifs_subset.txt.gz", package = "spatzie") motifs <- TFBSTools::readJASPARMatrix(motifs_file, matrixClass = "PFM") left <- GenomicRanges::GRanges( seqnames = c("chr1", "chr1", "chr1"), ranges = IRanges::IRanges(start = c(1, 15, 20), end = c(10, 35, 31))) right <- GenomicRanges::GRanges( seqnames = c("chr1", "chr2", "chr2"), ranges = IRanges::IRanges(start = c(17, 47, 41), end = c(28, 54, 53))) test_interactions <- GenomicInteractions::GenomicInteractions(left, right) # toy DNAStringSet to replace BSgenome object seqs <- c("chr1" = "CCACTAGCCACGCGTCACTGGTTAGCGTGATTGAAACTAAATCGTATGAAAATCC", "chr2" = "CTACAAACTAGGAATTTAGGCAAACCTGTGTTAAAATCTTAGCTCATTCATTAAT") toy_genome <- Biostrings::DNAStringSet(seqs, use.names = TRUE) res <- scan_motifs(test_interactions, motifs, toy_genome)
Looks for motifs which are significantly co-enriched from enhancer-promoter interaction data, derived from assays such as as HiC, ChIA-PET, etc. It can also look for differentially enriched motif pairs between to interaction experiments.
Jennifer Hammelman
Konstantin Krismer