Abstract
isomiRs package version: 1.10.1Lorena Pantano Harvard TH Chan School of Public Health, Boston, US
miRNAs are small RNA fragments (18-23 nt long) that influence gene expression during development and cell stability. Morin et al (Morin et al. 2008), discovered isomiRs first time after sequencing human stem cells.
IsomiRs are miRNAs that vary slightly in sequence, which result from variations in the cleavage site during miRNA biogenesis (5’-trimming and 3’-trimming variants), nucleotide additions to the 3’-end of the mature miRNA (3’-addition variants) and nucleotide modifications (substitution variants)(Martí et al. 2010).
There are many tools designed for isomiR detection, however the majority are web application where user can not control the analysis. The two main command tools for isomiRs mapping are SeqBuster and sRNAbench (Guillermo et al. 2014). isomiRs. package is designed to analyze the output of SeqBuster tool or any other tool after converting to the desire format.
If you use the package, please cite this paper (Pantano L 2010).
The input should be the output of SeqBuster-miraligner tool (*.mirna files). It is compatible with mirTOP tool as well, which parses BAM files with alignments against miRNA precursors.
For each sample the file should have the following format:
seq name freq mir start end mism add t5 t3
TGTAAACATCCTACACTCAGCT seq_100014_x23 23 hsa-miR-30b-5p 17 40 0 0 0 GT
TGTAAACATCCCTGACTGGAA seq_100019_x4 4 hsa-miR-30d-5p 6 26 13TC 0 0 g
TGTAAACATCCCTGACTGGAA seq_100019_x4 4 hsa-miR-30e-5p 17 37 12CT 0 0 g
CAAATTCGTATCTAGGGGATT seq_100049_x1 1 hsa-miR-10a-3p 63 81 0 TT 0 ata
TGACCTAGGAATTGACAGCCAGT seq_100060_x1 1 hsa-miR-192-5p 25 47 8GT 0 c agt
This is the standard output of SeqBuster-miraligner tool, but can be converted from any other tool having the mapping information on the precursors. Read more on miraligner manual.
This object will store all raw data from the input files and some processed information used for visualization and statistical analysis. It is a subclass of SummarizedExperiment
with colData
and counts
methods. Beside that, the object contains raw and normalized counts from miraligner allowing to update the summarization of miRNA expression.
The user can access the normalized count matrix with counts(object, norm=TRUE)
.
You can browse for the same miRNA or isomiRs in all samples with isoSelect
method.
## DataFrame with 6 rows and 14 columns
## cc1 cc2 cc3 cc4 cc5
## <numeric> <numeric> <numeric> <numeric> <numeric>
## hsa-let-7a-5p 235825 171354 149541 180654 168884
## hsa-let-7a-5p;iso_3p:T 8806 5478 5427 6249 5538
## hsa-let-7a-5p;iso_3p:gtt 1173 884 751 1010 963
## hsa-let-7a-5p;iso_3p:t 50920 52403 35525 53213 51499
## hsa-let-7a-5p;iso_3p:tt 6041 5467 3817 7491 6394
## hsa-let-7a-5p;iso_add:A 9616 4742 5437 6714 5846
## cc6 cc7 ct1 ct2 ct3
## <numeric> <numeric> <numeric> <numeric> <numeric>
## hsa-let-7a-5p 107430 153061 143030 163569 114028
## hsa-let-7a-5p;iso_3p:T 3734 5482 4825 6045 4626
## hsa-let-7a-5p;iso_3p:gtt 580 863 927 756 495
## hsa-let-7a-5p;iso_3p:t 30659 45492 42878 35795 24993
## hsa-let-7a-5p;iso_3p:tt 3984 5940 7107 4828 2858
## hsa-let-7a-5p;iso_add:A 3799 5532 5214 7326 4790
## ct4 ct5 ct6 ct7
## <numeric> <numeric> <numeric> <numeric>
## hsa-let-7a-5p 123454 133092 158909 140272
## hsa-let-7a-5p;iso_3p:T 4181 4505 5134 4670
## hsa-let-7a-5p;iso_3p:gtt 719 682 1048 860
## hsa-let-7a-5p;iso_3p:t 34512 36973 47873 39585
## hsa-let-7a-5p;iso_3p:tt 5087 5652 6589 5447
## hsa-let-7a-5p;iso_add:A 4430 5010 6034 5081
metadata(mirData)
contains two lists: rawList
is a list with same length than number of samples and stores the input files for each sample; isoList
is a list with same length than number of samples and stores information for each isomiR type summarizing the different changes for the different isomiRs (trimming at 3’, trimming a 5’, addition and substitution). For instance, you can get the data stored in isoList
for sample 1 and 5’ changes with this code metadata(ids)[["isoList"]][[1]]["t5sum"]
.
IsomiR names follows this structure:
iso
if the sequence has variations. t5 tag: indicates variations at 5’ position. The naming contains two words: direction - nucleotides
, where direction can be UPPER CASE NT (changes upstream of the 5’ reference position) or LOWER CASE NT (changes downstream of the 5’ reference position). 0
indicates no variation, meaning the 5’ position is the same than the reference. After direction
, it follows the nucleotide/s that are added (for upstream changes) or deleted (for downstream changes). t3 tag: indicates variations at 3’ position. The naming contains two words: direction - nucleotides
, where direction can be LOWER CASE NT (upstream of the 3’ reference position) or UPPER CASE NT (downstream of the 3’ reference position). 0
indicates no variation, meaning the 3’ position is the same than the reference. After direction
, it follows the nucleotide/s that are added (for downstream changes) or deleted (for upstream chanes). ad tag: indicates nucleotides additions at 3’ position. The naming contains two words: direction - nucleotides
, where direction is UPPER CASE NT (upstream of the 5’ reference position). 0
indicates no variation, meaning the 3’ position has no additions. After direction
, it follows the nucleotide/s that are added. mm tag: indicates nucleotides substitutions along the sequences. The naming contains three words: position-nucleotideATsequence-nucleotideATreference
. *seed tag: same than mm
tag, but only if the change happens between nucleotide 2 and 8.In general nucleotides in UPPER case mean insertions respect to the reference sequence, and nucleotides in LOWER case mean deletions respect to the reference sequence.
We are going to use a small RNAseq data from human brain samples (???) to give some basic examples of isomiRs analyses.
In this data set we will find two groups:
pc: 7 control individuals pt: 7 patients with Parkinson’s Disease in early stage.
The function IsomirDataSeqFromFiles
needs a vector with the paths for each file and a data frame with the design experiment similar to the one used for a mRNA differential expression analysis. Row names of the data frame should be the names for each sample in the same order than the list of files.
You can plot isomiRs expression with isoPlot
. In this figure you will see how abundant is each type of isomiRs at different positions considering the total abundance and the total number of sequences. The type
parameter controls what type of isomiRs to show. It can be trimming (iso5 and iso3), addition (add) or substitution (subs) changes.
isoCounts
gets the count matrix that can be used for many different downstream analyses changing the way isomiRs are collapsed. The following command will merge all isomiRs into one feature: the reference miRNA.
## cc1 cc2 cc3 cc4 cc5 cc6 cc7 ct1
## hsa-let-7a-2-3p 5 13 4 9 11 6 7 0
## hsa-let-7a-3p 767 707 630 609 731 389 681 258
## hsa-let-7a-5p 317730 244832 203888 260931 244784 153049 220563 208511
## hsa-let-7b-3p 1037 1949 679 1385 1884 697 1499 338
## hsa-let-7b-5p 69671 64702 42347 65322 60767 34975 56875 22215
## hsa-let-7c-3p 36 16 25 40 26 23 35 3
## ct2 ct3 ct4 ct5 ct6 ct7
## hsa-let-7a-2-3p 4 2 12 12 8 3
## hsa-let-7a-3p 496 343 635 622 452 519
## hsa-let-7a-5p 222066 154246 175523 189733 230690 199923
## hsa-let-7b-3p 931 494 1149 1194 1431 988
## hsa-let-7b-5p 45069 28312 41214 43058 61258 46600
## hsa-let-7c-3p 30 28 26 24 22 17
The normalization uses rlog
from DESeq2 package and allows quick integration to another analyses like heatmap, clustering or PCA.
library(pheatmap)
ids = isoNorm(ids, formula = ~ condition)
pheatmap(counts(ids, norm=TRUE)[1:100,],
annotation_col = data.frame(colData(ids)[,1,drop=FALSE]),
show_rownames = FALSE, scale="row")
The isoDE
uses functions from DESeq2 package. This function has parameters to create a matrix using only the reference miRNAs, all isomiRs, or some of them. This matrix and the design matrix are the inputs for DESeq2. The output will be a DESeqDataSet object, allowing to generate any plot or table explained in DESeq2 package vignette.
## log2 fold change (MLE): condition ct vs cc
## Wald test p-value: condition ct vs cc
## DataFrame with 6 rows and 6 columns
## baseMean log2FoldChange lfcSE
## <numeric> <numeric> <numeric>
## hsa-let-7a-2-3p 6.60005242130612 -0.391808034451504 0.609492802121228
## hsa-let-7a-3p 546.805469613238 -0.225685812620638 0.224797224087315
## hsa-let-7a-5p 220625.551423747 0.113470209115605 0.270505352954979
## hsa-let-7b-3p 1096.78693283053 -0.391274509777866 0.329460850713291
## hsa-let-7b-5p 47415.9668329051 -0.235214733416519 0.226064546587921
## hsa-let-7c-3p 23.1177591558883 -0.174438078886159 0.302732297642984
## stat pvalue padj
## <numeric> <numeric> <numeric>
## hsa-let-7a-2-3p -0.642842758910176 0.520326134785159 0.989527585749852
## hsa-let-7a-3p -1.00395284477791 0.315401343211605 0.989527585749852
## hsa-let-7a-5p 0.419474911960393 0.674869086101299 0.989527585749852
## hsa-let-7b-3p -1.18762065031626 0.234982898335359 0.989527585749852
## hsa-let-7b-5p -1.04047599221862 0.298118812264401 0.989527585749852
## hsa-let-7c-3p -0.576212317761602 0.564471680375386 0.989527585749852
You can differentiate between reference sequences and isomiRs at 5’ end with this command:
## row baseMean log2FoldChange lfcSE
## 1 hsa-let-7a-2-3p 2.4059793 -0.8425278 1.0440442
## 2 hsa-let-7a-2-3p;iso_5p:TAA 0.1285419 1.1051761 3.0845751
## 3 hsa-let-7a-2-3p;ref 4.0407844 -0.2219555 0.6581561
## 4 hsa-let-7a-3p 376.7216287 -0.1663188 0.2027009
## 5 hsa-let-7a-3p;iso_5p:A 0.6517381 1.0223121 1.7714731
## 6 hsa-let-7a-3p;iso_5p:c 84.6025648 -0.3165339 0.2538101
## stat pvalue padj
## 1 -0.8069848 0.4196753 0.9865695
## 2 0.3582912 0.7201254 0.9865695
## 3 -0.3372384 0.7359372 0.9865695
## 4 -0.8205135 0.4119234 0.9865695
## 5 0.5770972 0.5638738 0.9865695
## 6 -1.2471289 0.2123502 0.9865695
Alternative, for more complicated cases or if you want to control more the differential expression analysis paramters you can use directly DESeq2 package feeding it with the output of counts(ids)
and colData(ids)
like this:
Partial Least Squares Discriminant Analysis (PLS-DA) is a technique specifically appropriate for analysis of high dimensionality data sets and multicollineality (Pérez-Enciso and Tenenhaus 2003). PLS-DA is a supervised method (i.e. makes use of class labels) with the aim to provide a dimension reduction strategy in a situation where we want to relate a binary response variable (in our case young or old status) to a set of predictor variables. Dimensionality reduction procedure is based on orthogonal transformations of the original variables (isomiRs) into a set of linearly uncorrelated latent variables (usually termed as components) such that maximizes the separation between the different classes in the first few components (Xia and Wishart 2011). We used sum of squares captured by the model (R2) as a goodness of fit measure. We implemented this method using the DiscriMiner
into isoPLSDA
function. The output p-value of this function will tell about the statistical significant of the group separation using miRNA expression data. Moreover, the function isoPLSDAplot
helps to visualize the results. It will plot the samples using the significant components (t1, t2, t3 …) from the PLS-DA analysis and the samples distribution along the components.
ids = isoCounts(ids, iso5=TRUE, minc=10, mins=6)
ids = isoNorm(ids, formula = ~ condition)
pls.ids = isoPLSDA(ids, "condition", nperm = 2)
df = isoPLSDAplot(pls.ids)
The analysis can be done again using only the most important discriminant isomiRS from the PLS-DA models based on the analysis. We used Variable Importance for the Projection (VIP) criterion to select the most important features, since takes into account the contribution of a specific predictor for both the explained variability on the response and the explained variability on the predictors.
The package offers a correlation analysis of miRNA and gene expression data. Having two SummarizedExperiments objects with their expression, the target prediction for each miRNA, the function isoNetwork
and isoPlotNetwork
can generate a summarized figure showing the relationship between expression profile, miRNA repression and enrichment analysis:
# library(org.Mm.eg.db)
# library(clusterProfiler)
data(isoExample)
# ego <- enrichGO(row.names(assay(gene_ex_rse, "norm")),
# org.Mm.eg.db, "ENSEMBL", ont = "BP")
data = isoNetwork(mirna_ex_rse, gene_ex_rse, target = ma_ex,
enrich = ego, summarize = "group")
## Dimmension of cor matrix: 20 25
As an option, org
can be org.Mm.eg.db
and genename
can be ENSEMBL
and it will run enrcihGO
internally.
To create the ma_ex matrix, the user can use findTargets
:
mirna_ma <- data.frame(gene = names(gene_ex_rse)[1:20],
mir = names(mirna_ex_rse))
ma_ex <- findTargets(mirna_ex_rse, gene_ex_rse, mirna_ma)
## Dimmension of cor matrix: 20 20
## mmu-miR-101a-3p mmu-miR-107-3p mmu-miR-142a-3p
## ENSMUSG00000004530 0 0 0
## ENSMUSG00000005397 0 0 0
## ENSMUSG00000006386 0 0 0
## ENSMUSG00000016494 0 0 0
## ENSMUSG00000019997 0 0 0
## ENSMUSG00000021822 0 0 0
## mmu-miR-144-3p
## ENSMUSG00000004530 0
## ENSMUSG00000005397 0
## ENSMUSG00000006386 0
## ENSMUSG00000016494 0
## ENSMUSG00000019997 0
## ENSMUSG00000021822 0
And to get the mirna_ma
data.frame with the miRNA-target information, the user can use mirna2targetscan
function:
## miRFamily Seedmatch PCT gene targetscan
## 1 miR-34ac/34bc-5p/449abc/449c-5p 8mer 0.29 79719 hsa-miR-34c-5p
## 2 miR-34ac/34bc-5p/449abc/449c-5p 8mer 0.76 215 hsa-miR-34c-5p
## 3 miR-34ac/34bc-5p/449abc/449c-5p 7mer-1a 0.58 3983 hsa-miR-34c-5p
## 4 miR-34ac/34bc-5p/449abc/449c-5p 8mer 0.73 29 hsa-miR-34c-5p
## 5 miR-34ac/34bc-5p/449abc/449c-5p 8mer 0.67 64746 hsa-miR-34c-5p
## 6 miR-34ac/34bc-5p/449abc/449c-5p 8mer 0.75 40 hsa-miR-34c-5p
## mir ext num letter
## 1 hsa-miR-34c-5p hsa-miR-34c hsa-miR-34c hsa-miR-34
## 2 hsa-miR-34c-5p hsa-miR-34c hsa-miR-34c hsa-miR-34
## 3 hsa-miR-34c-5p hsa-miR-34c hsa-miR-34c hsa-miR-34
## 4 hsa-miR-34c-5p hsa-miR-34c hsa-miR-34c hsa-miR-34
## 5 hsa-miR-34c-5p hsa-miR-34c hsa-miR-34c hsa-miR-34
## 6 hsa-miR-34c-5p hsa-miR-34c hsa-miR-34c hsa-miR-34
Here is the output of sessionInfo
on the system on which this document was compiled:
## R version 3.5.2 (2018-12-20)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 16.04.5 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.8-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.8-bioc/R/lib/libRlapack.so
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_US.UTF-8 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] parallel stats4 stats graphics grDevices utils datasets
## [8] methods base
##
## other attached packages:
## [1] targetscan.Hs.eg.db_0.6.1 AnnotationDbi_1.44.0
## [3] DESeq2_1.22.2 pheatmap_1.0.12
## [5] bindrcpp_0.2.2 isomiRs_1.10.1
## [7] SummarizedExperiment_1.12.0 DelayedArray_0.8.0
## [9] BiocParallel_1.16.5 matrixStats_0.54.0
## [11] Biobase_2.42.0 GenomicRanges_1.34.0
## [13] GenomeInfoDb_1.18.1 IRanges_2.16.0
## [15] S4Vectors_0.20.1 BiocGenerics_0.28.0
## [17] DiscriMiner_0.1-29 BiocStyle_2.10.0
##
## loaded via a namespace (and not attached):
## [1] assertive.base_0.0-7 colorspace_1.3-2
## [3] rjson_0.2.20 circlize_0.4.5
## [5] htmlTable_1.13 XVector_0.22.0
## [7] GlobalOptions_0.1.0 base64enc_0.1-3
## [9] ggdendro_0.1-20 rstudioapi_0.8
## [11] assertive.sets_0.0-3 ggrepel_0.8.0
## [13] bit64_0.9-7 codetools_0.2-16
## [15] splines_3.5.2 logging_0.8-104
## [17] mnormt_1.5-5 geneplotter_1.60.0
## [19] knitr_1.21 Formula_1.2-3
## [21] Nozzle.R1_1.1-1 broom_0.5.1
## [23] annotate_1.60.0 cluster_2.0.7-1
## [25] readr_1.3.1 BiocManager_1.30.4
## [27] compiler_3.5.2 backports_1.1.3
## [29] assertthat_0.2.0 Matrix_1.2-15
## [31] lazyeval_0.2.1 limma_3.38.3
## [33] lasso2_1.2-20 acepack_1.4.1
## [35] htmltools_0.3.6 tools_3.5.2
## [37] gtable_0.2.0 glue_1.3.0
## [39] GenomeInfoDbData_1.2.0 reshape2_1.4.3
## [41] dplyr_0.7.8 Rcpp_1.0.0
## [43] gdata_2.18.0 nlme_3.1-137
## [45] psych_1.8.10 xfun_0.4
## [47] stringr_1.3.1 DEGreport_1.18.1
## [49] gtools_3.8.1 XML_3.98-1.16
## [51] edgeR_3.24.3 zlibbioc_1.28.0
## [53] MASS_7.3-51.1 scales_1.0.0
## [55] hms_0.4.2 RColorBrewer_1.1-2
## [57] ComplexHeatmap_1.20.0 yaml_2.2.0
## [59] memoise_1.1.0 gridExtra_2.3
## [61] ggplot2_3.1.0 rpart_4.1-13
## [63] reshape_0.8.8 latticeExtra_0.6-28
## [65] stringi_1.2.4 RSQLite_2.1.1
## [67] genefilter_1.64.0 checkmate_1.8.5
## [69] caTools_1.17.1.1 shape_1.4.4
## [71] rlang_0.3.0.1 pkgconfig_2.0.2
## [73] bitops_1.0-6 evaluate_0.12
## [75] lattice_0.20-38 purrr_0.2.5
## [77] bindr_0.1.1 labeling_0.3
## [79] htmlwidgets_1.3 cowplot_0.9.3
## [81] bit_1.1-14 tidyselect_0.2.5
## [83] GGally_1.4.0 plyr_1.8.4
## [85] magrittr_1.5 R6_2.3.0
## [87] gplots_3.0.1 generics_0.0.2
## [89] Hmisc_4.1-1 DBI_1.0.0
## [91] pillar_1.3.1 foreign_0.8-71
## [93] survival_2.43-3 RCurl_1.95-4.11
## [95] nnet_7.3-12 tibble_2.0.0
## [97] crayon_1.3.4 KernSmooth_2.23-15
## [99] rmarkdown_1.11 GetoptLong_0.1.7
## [101] locfit_1.5-9.1 grid_3.5.2
## [103] data.table_1.11.8 blob_1.1.1
## [105] ConsensusClusterPlus_1.46.0 digest_0.6.18
## [107] xtable_1.8-3 tidyr_0.8.2
## [109] munsell_0.5.0
Guillermo, Barturen, Rueda Antonio, Hamberg Maarten, Alganza Angel, Lebron Ricardo, Kotsyfakis Michalis, Shi BuJun, KoppersLalic Danijela, and Hackenberg Michael. 2014. “sRNAbench: profiling of small RNAs and its sequence variants in single or multi-species high-throughput experiments.” Methods in Next Generation Sequencing 1 (1):2084–7173. https://doi.org/10.2478/mngs-2014-0001.
Martí, Eulàlia, Lorena Pantano, Mónica Bañez-Coronel, Franc Llorens, Elena Miñones-Moyano, Sílvia Porta, Lauro Sumoy, Isidre Ferrer, and Xavier Estivill. 2010. “A myriad of miRNA variants in control and Huntington’s disease brain regions detected by massively parallel sequencing.” Nucleic Acids Res. 38:7219–35. https://doi.org/10.1093/nar/gkq575.
Morin, R. D., M. D. O’Connor, M. Griffith, F. Kuchenbauer, A. Delaney, A.-L. Prabhu, Y. Zhao, et al. 2008. “Application of massively parallel sequencing to microRNA profiling and discovery in human embryonic stem cells.” Genome Res. 18:610–21. https://doi.org/10.1101/gr.7179508.
Pantano L, Estivil X, Marti E. 2010. “SeqBuster.” Nucleic Acids Res. 38:e34. https://doi.org/10.1093/nar/gkp1127.
Pérez-Enciso, Miguel, and Michel Tenenhaus. 2003. “Prediction of clinical outcome with microarray data: a partial least squares discriminant analysis (PLS-DA) approach.” Human Genetics 112:581–92. https://doi.org/10.1007/s00439-003-0921-9.
Xia, Jianguo, and David S Wishart. 2011. “Web-based inference of biological patterns, functions and pathways from metabolomic data using MetaboAnalyst.” Nature Protocols 6:743–60. https://doi.org/10.1038/nprot.2011.319.