scanMiR 1.2.0
## Warning in knitr::include_graphics(system.file("docs", "sticker.svg", package =
## "scanMiR")): It is highly recommended to use relative paths for images. You had
## absolute paths: "/tmp/RtmpMATF3z/Rinste54c68cda8d0/scanMiR/docs/sticker.svg"
scanMiR can be used to identify potential binding sites given a
set of miRNAs and a set of transcripts. Furthermore, it determines the type of
binding site and, given a KdModel
object, the predicted affinity of the site.
The main function used for determining matches of miRNAs in a given set of
sequences is findSeedMatches
. It accepts a set of DNA sequences either as a
character vector or as a DNAStringSet.
The miRNAs can be provided either as a character vector of seeds/miRNA
sequences or as a KdModelList
.
The seed must be given in the form of a (RNA or DNA) sequence of length 7 or 8 (the 8th nucleotide being the final ‘A’ - it is added if only 7 are given). Note that the seed should be given as it would appear in a match in the target sequence (i.e. the reverse complement of how it appears in the miRNA).
library(scanMiR)
# seed sequence of hsa-miR-155-5p
seed <- "AGCAUUAA"
# load a sample transcript
data("SampleTranscript")
# run scan
matches <- findSeedMatches(SampleTranscript, seed, verbose = FALSE)
matches
## GRanges object with 2 ranges and 1 metadata column:
## seqnames ranges strand | type
## <Rle> <IRanges> <Rle> | <factor>
## [1] seq1 491-498 * | 8mer
## [2] seq1 692-699 * | 7mer-m8
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
By default, a GRanges
object is returned. Apart from the position of the matches, it provides
information on the type of the putative binding site corresponding to the match.
Setting ret = "data.frame"
returns the same information as a data.frame.
Alternatively, we can provide the full miRNA sequence, which results in additional information on supplementary 3’ pairing in the form of an aggregated score (see Section 1.3.2 for further details).
# full sequence of the mature miR-155-5p transcript
miRNA <- "UUAAUGCUAAUCGUGAUAGGGGUU"
# run scan
matches <- findSeedMatches(SampleTranscript, miRNA, verbose = FALSE)
## Warning in XStringSet(pattern_seqtype, pattern): metadata columns on input
## DNAStringSet object were dropped
matches
## GRanges object with 2 ranges and 3 metadata columns:
## seqnames ranges strand | type p3.score note
## <Rle> <IRanges> <Rle> | <factor> <integer> <Rle>
## [1] seq1 491-498 * | 8mer 12 TDMD?
## [2] seq1 692-699 * | 7mer-m8 0 -
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
We can take a closer look at the alignment of the first match:
viewTargetAlignment(matches[1], miRNA, SampleTranscript)
## Warning in XStringSet(pattern_seqtype, pattern): metadata columns on input
## DNAStringSet object were dropped
##
## miRNA 3'-UUGGGGAUAGUGC-UAA-UCGUAAUU-5'
## |||||||||||| ||||||||
## target 5'-...NGAGUAACGACCCCUAUCACGUCCGCAGCAUUAAAU...-3'
Apart from the direct seed match (right), this representation also reveals the extensive supplementary 3’ pairing (left).
Finally, we can provide the miRNA in the form of a KdModel
(see the
vignette on KdModels for further information). In this case
findSeedMatches
also returns the predicted affinity value for each match. The
log_kd
column contains log(Kd) values multiplied by 1000, where Kd is the
predicted dissociation constant of miRNA:mRNA binding for the putative binding
site.
# load sample KdModel
data("SampleKdModel")
# run scan
matches <- findSeedMatches(SampleTranscript, SampleKdModel, verbose = FALSE)
## Warning in XStringSet(pattern_seqtype, pattern): metadata columns on input
## DNAStringSet object were dropped
matches
## GRanges object with 2 ranges and 4 metadata columns:
## seqnames ranges strand | type log_kd p3.score note
## <Rle> <IRanges> <Rle> | <factor> <integer> <integer> <Rle>
## [1] seq1 491-498 * | 8mer -4868 12 TDMD?
## [2] seq1 692-699 * | 7mer-m8 -3702 0 -
## -------
## seqinfo: 1 sequence from an unspecified genome; no seqlengths
Running findSeedMatches
using a Kdmodel
also returns matches that correspond
to non-canonical binding sites. These are typically of low affinity, but may
affect repression if several of them are found on the same transcript. The scan
can be restricted to canonical sites using the option onlyCanonical = TRUE
.
KdModel
collections corresponding to all human, mouse and rat mirbase miRNAs
can be obtained through the
scanMiRData package.
If the transcript sequences are provided as a DNAStringSet
, one may specify
the length of the open reading frame region of the transcripts as a metadata
column in order to distinguish between matches in the ORF and 3’UTR regions.
library(Biostrings)
# generate set of random sequences
seqs <- DNAStringSet(getRandomSeq(length = 1000, n = 10))
# add vector of ORF lengths
mcols(seqs)$ORF.length <- sample(500:800, length(seqs))
# run scan
matches2 <- findSeedMatches(seqs, SampleKdModel, verbose = FALSE)
## Warning in XStringSet(pattern_seqtype, pattern): metadata columns on input
## DNAStringSet object were dropped
head(matches2)
## GRanges object with 4 ranges and 5 metadata columns:
## seqnames ranges strand | type log_kd p3.score note ORF
## <Rle> <IRanges> <Rle> | <factor> <integer> <integer> <Rle> <Rle>
## [1] seq3 838-845 * | non-canonical -1002 0 - FALSE
## [2] seq10 745-752 * | 6mer -3336 0 - TRUE
## [3] seq10 763-770 * | 6mer -1933 4 - TRUE
## [4] seq10 452-459 * | non-canonical -1545 6 - TRUE
## -------
## seqinfo: 10 sequences from an unspecified genome; no seqlengths
Upon binding the seed regions, further complementary pairing of the target to the 3’
region of the miRNA can increase affinity and further stabilize the binding
(Schirle, Sheu-Gruttadauria and MacRae, 2014).
Upon finding a seed match, scanMiR
performs a local alignment on the upstream
region to identify such complementary 3’ binding. This is internally done by the
get3pAlignment()
function, the arguments to which (e.g. the maximum size of
the gap between the seed binding and the complementary binding) can be passed
via the findSeedMatches
argument p3.params
. By default, when running
findSeedMatches
a 3’ score is reported in the matches, which roughly
corresponds to the number of consecutive matching nucleotides (adjusting for
small gaps and T/G substitutions) within the constraints (see
?get3pAlignment
for more detail). More information (such as the size of the
miRNA and target loops between the two complementary regions) can be reported by
setting findSeedMatches(..., p3.extra=TRUE)
. In addition, the pairing can be
visualized with viewTargetAlignment
:
viewTargetAlignment(matches[1], SampleKdModel, SampleTranscript)
## Warning in XStringSet(pattern_seqtype, pattern): metadata columns on input
## DNAStringSet object were dropped
##
## miRNA 3'-UUGGGGAUAGUGC-UAA-UCGUAAUU-5'
## |||||||||||| ||||||||
## target 5'-...NGAGUAACGACCCCUAUCACGUCCGCAGCAUUAAAU...-3'
Some forms of 3’ bindings can however lead to drastic functional consequences.
For example, sufficient final complementary at the 3’ end of the miRNA can lead
to Target-Directed miRNA Degradation (TDMD,
Sheu-Gruttadauria, Pawlica et al., 2019).
findSeedMatches
will also flag such putative sites in the notes
column of
the matches. Finally, while circular RNAs can act as miRNA sponges, some miRNA
bindings can slice their circular structure
Hansen et al., 2011 and free their
cargo. findSeedMatches
will also flag such sites in the notes
column.
The shadow
argument can be used to take into account the observation that
sites within the first ~15 nucleotides of the 3’UTR show poor efficiency
(Grimson et al. 2007).
findSeedMatches
will treat matches within the first shadow
positions of the
UTR in the same way as matches in the ORF region. If no information on ORF
lengths is provided, it will simply ignore the first shadow
positions. The
default setting is shadow = 0
.
The parameter minDist
can be used to specify the minimum distance between
matches of the same miRNA (default 7). If there are multiple matches within
minDist
, only the highest affinity match will be considered.
With ret = "aggregated"
one obtains a data.frame that contains the predicted
repression for each sequence-seed-pair aggregated over all matches along with
information about the types and number of matches. Parameters for aggregation
can be specified using agg.params
. For further details, see Section
2.
scanMiR implements aggregation of miRNA sites based on the biochemical model from McGeary, Lin et al. (2019). This model first predicts the occupancy of AGO-miRNA complexes at each potential binding site as a function of the measured or estimated dissociation constants (Kds). It then assumes an additive effect of the miRNA on the basal decay rate of the transcript that is proportional to this occupancy.
The key parameters of this model are:
a
: the relative concentration of unbound AGO-miRNA complexesb
: the factor that multiplies with the occupancy and is added to the basal
decay rate (can be interpreted as the additional repression caused by a single
bound AGO)c
: the penalty factor for sites that are found within the ORF regionMore specifically, the occupancy of a mRNA \(m\) by miRNA \(g\), with \(p\) matches in the ORF region and \(q\) matches in the 3’UTR region, is given by the following equation: \[ \begin{equation} N_{m,g} = \sum_{i=1}^{p}\left(\frac{a_g}{a_g + c_{\text{ORF}} K_{d,i}^{\text{ORF}}}\right) + \sum_{j=1}^{q}\left(\frac{a_g}{a_g + K_{d,j}^{\text{3'UTR}}}\right) \end{equation} \]
The corresponding background occupancy is estimated by substituting the average
affinity of nonspecifically bound sites (i.e. \(K_d = 1.0\)):
\[
\begin{equation}
N_{m,g,\text{background}} =
\sum_{i=1}^{p}\left(\frac{a_g}{a_g + c_{\text{ORF}}}\right) +
\sum_{j=1}^{q}\left(\frac{a_g}{a_g + 1}\right)
\end{equation}
\]
In addition to this original model, scanMiR
includes a coefficient e
which
adjusts the Kd values based on the supplementary 3’ alignment:
\[ \begin{equation} N_{m,g} = \sum_{i=1}^{p}\left(\frac{a_g}{a_g + e_{i}c_{\text{ORF}} K_{d,i}^{\text{ORF}}}\right) + \sum_{j=1}^{q}\left(\frac{a_g}{a_g + e_{j}K_{d,j}^{\text{3'UTR}}}\right) \end{equation} \]
with \(e_i = \exp(\text{p3}\cdot\text{p3.score}_i)\). p3
is a global parameter,
and \(p3.score_i\) is the 3’ alignment score (roughly corresponding to the number
of matched nucleotides, by default capped to 8 and set to 0 if below 3). Of
note, the default value of p3
is very small, leading to a very mild effect.
The importance of complementary binding seems to depend on the miRNA, and at the
moment there is no easy way to predict this from the miRNA sequence. Our
conservative estimate might not attribute sufficient importance to this factor
for some miRNAs.
The repression is then obtained as the log fold change of the two occupancies: \[ \text{repression} = \log(1+bN_{m,g,\text{background}}) - \log(1+bN_{m,g}) \]
Because UTR and ORF lengths have been reported to influence the efficacy of
repression, scanMiR
also includes an additional modifier to terms handling
these effects:
\[
\text{repression}_{\text{adj}} = \text{repression}\cdot
(1+f\cdot\text{UTR.length}+h\cdot\text{ORF.length})
\]
While b
, c
, p3
, f
and h
are considered global parameters (i.e. the
same for different miRNAs and transcripts and also across experimental
contexts), a
is expected to be different for each miRNA in a
given experimental condition. However, as shown by
McGeary, Lin et al. (2019), the
model performance is robust to changes in a
over several orders of magnitude.
Aggregation for all miRNA-transcript pairs for a given data set is therefore
usually based on a single a
value.
Given a GRanges
or data.frame of matches as returned by findSeedMatches
,
aggregation can be performed by the function aggregateMatches
:
agg_matches <- aggregateMatches(matches2)
head(agg_matches)
## transcript repression 8mer 7mer 6mer non-canonical ORF.canonical
## 1 seq10 -0.12609016 0 0 0 0 2
## 2 seq3 -0.03232475 0 0 0 1 0
## ORF.non-canonical
## 1 1
## 2 0
This returns a data.frame with predicted repression values for each
miRNA-transcript pair along with a count table of the different site types. If
matches
does not contain a log_kd
column, only the count table will be
returned.
scanMiR uses the following default parameter values for aggregation that have been determined by fitting and validating the model using several experimental data sets:
unlist(scanMiR:::.defaultAggParams())
## a b c p3 coef_utr coef_orf
## 0.007726 0.573500 0.181000 0.051000 -0.171060 -0.215460
## keepSiteInfo
## 1.000000
Where coef_utr
and coef_orf
respectively correspond to the f
and h
in
the above formula. To disable these features, they can simply be set to
zero. keepSiteInfo
lets you choose whether the site count table should be
returned. The parameters can be passed directly to aggregateMatches
, or passed
to findSeedMatch
when doing aggregation on-the-fly using the agg.params
argument.
To deal with large amounts of sequences and/or seeds, both findSeedMatches
and aggregateMatches
support multithreading using the
BiocParallel package. This can be activated by passing
BP = MulticoreParam(ncores)
.
Depending on the system and the size of the scan (i.e. when including all
non-canonical sites), mutlithreading can potentially take a large amount of
memory. To avoid memory issues, the number of seeds processed simultaneously by
findSeedMatches
can be restricted using the n_seeds
parameter.
Alternatively, scan results can be saved to temporary files using the
useTmpFiles
argument (see ?findSeedMatches
for more detail).
Note that in addition to the multithreading specified in its arguments,
aggregateMatches
uses the data.table package, which is often
set to use multi-threading by default (see ?data.table::setDTthreads
for more
information). This can leave to CPU usage higher than specified through the
BP
argument of aggregateMatches
.
Binding sites for all miRNAs on all transcripts, especially when including
non-canonical sites, can easily amount to prohibitive amounts of memory. The
companion scanMiRApp package includes
a class implementing fast indexed access to on-disk GenomicRanges and
data.frames. The package additionally contains wrapper (e.g. for performing full
transcriptome scans) for common species and for detecting enriched miRNA-target
pairs, as well as a shiny interface to scanMiR
.
## R version 4.2.0 RC (2022-04-19 r82224)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.4 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.15-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.15-bioc/R/lib/libRlapack.so
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_GB LC_COLLATE=C
## [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
## [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
## [9] LC_ADDRESS=C LC_TELEPHONE=C
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## attached base packages:
## [1] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] Biostrings_2.64.0 GenomeInfoDb_1.32.0 XVector_0.36.0
## [4] IRanges_2.30.0 S4Vectors_0.34.0 BiocGenerics_0.42.0
## [7] scanMiR_1.2.0 BiocStyle_2.24.0
##
## loaded via a namespace (and not attached):
## [1] tidyselect_1.1.2 xfun_0.30 bslib_0.3.1
## [4] purrr_0.3.4 colorspace_2.0-3 vctrs_0.4.1
## [7] generics_0.1.2 htmltools_0.5.2 yaml_2.3.5
## [10] utf8_1.2.2 rlang_1.0.2 jquerylib_0.1.4
## [13] pillar_1.7.0 glue_1.6.2 DBI_1.1.2
## [16] BiocParallel_1.30.0 GenomeInfoDbData_1.2.8 lifecycle_1.0.1
## [19] ggseqlogo_0.1 stringr_1.4.0 zlibbioc_1.42.0
## [22] munsell_0.5.0 gtable_0.3.0 evaluate_0.15
## [25] labeling_0.4.2 knitr_1.38 fastmap_1.1.0
## [28] parallel_4.2.0 fansi_1.0.3 highr_0.9
## [31] Rcpp_1.0.8.3 scales_1.2.0 BiocManager_1.30.17
## [34] magick_2.7.3 jsonlite_1.8.0 farver_2.1.0
## [37] ggplot2_3.3.5 digest_0.6.29 stringi_1.7.6
## [40] bookdown_0.26 dplyr_1.0.8 cowplot_1.1.1
## [43] GenomicRanges_1.48.0 grid_4.2.0 cli_3.3.0
## [46] tools_4.2.0 bitops_1.0-7 magrittr_2.0.3
## [49] sass_0.4.1 RCurl_1.98-1.6 tibble_3.1.6
## [52] crayon_1.5.1 pkgconfig_2.0.3 ellipsis_0.3.2
## [55] data.table_1.14.2 assertthat_0.2.1 rmarkdown_2.14
## [58] R6_2.5.1 compiler_4.2.0