mitch is an R package for multi-contrast enrichment analysis. At it’s heart, it uses a rank-MANOVA based statistical approach to detect sets of genes that exhibit enrichment in the multidimensional space as compared to the background. The rank-MANOVA concept dates to work by Cox and Mann (https://doi.org/10.1186/1471-2105-13-S16-S12). mitch is useful for pathway analysis of profiling studies with one, two or more contrasts, or in studies with multiple omics profiling, for example proteomic, transcriptomic, epigenomic analysis of the same samples. mitch is perfectly suited for pathway level differential analysis of scRNA-seq data.
The main strengths of mitch are that it can import datasets easily from many upstream tools and has advanced plotting features to visualise these enrichments. mitch consists of five functions. A typical mitch workflow would consist of:
1: Import gene sets with gmt_import()
2: Import profiling data with mitch_import()
3: Calculate enrichments with mitch_calc()
4: And generate plots and reports with mitch_plots()
, mitch_report()
and networkplot()
.
if(!requireNamespace("BiocManager", quietly = TRUE))
install.packages("BiocManager")
::install("mitch") BiocManager
library("mitch")
mitch has a function to import GMT files to R lists (adapted from Yu et al, 2012 in the clusterProfiler package). For example we can grab some gene sets from Reactome.org:
download.file("https://reactome.org/download/current/ReactomePathways.gmt.zip",
destfile="ReactomePathways.gmt.zip")
unzip("ReactomePathways.gmt.zip")
gmt_import("ReactomePathways.gmt") genesets <-
In this cut down example we will be using a sample of 200 Reactome gene sets:
data(genesetsExample)
head(genesetsExample,3)
## $`2-LTR circle formation`
## [1] "Reactome Pathway" "BANF1" "HMGA1" "LIG4"
## [5] "PSIP1" "XRCC4" "XRCC5" "XRCC6"
## [9] "gag" "gag-pol" "rev" "vif"
## [13] "vpr" "vpu"
##
## $`5-Phosphoribose 1-diphosphate biosynthesis`
## [1] "Reactome Pathway" "PRPS1" "PRPS1L1" "PRPS2"
##
## $`A tetrasaccharide linker sequence is required for GAG synthesis`
## [1] "Reactome Pathway" "AGRN" "B3GALT6" "B3GAT1"
## [5] "B3GAT2" "B3GAT3" "B4GALT7" "BCAN"
## [9] "BGN" "CSPG4" "CSPG5" "DCN"
## [13] "GPC1" "GPC2" "GPC3" "GPC4"
## [17] "GPC5" "GPC6" "HSPG2" "NCAN"
## [21] "SDC1" "SDC2" "SDC3" "SDC4"
## [25] "VCAN" "XYLT1" "XYLT2"
mitch accepts pre-ranked data supplied by the user, but also has a function called mitch_import
for importing tables generated by limma, edgeR, DESeq2, ABSSeq, Sleuth, Seurat, Muscat and several other upstream tools. By default, only the genes that are detected in all contrasts are included, but this behaviour can be modified for sparse data setting joinType=full
. The below example imports two edgeR tables called “rna” and “k9a” Where gene identifiers are present as row names. Note that if there is more than one profile being imported, they need to be part of a list.
data(rna,k9a)
list("rna"=rna,"k9a"=k9a)
x <- mitch_import(x,"edgeR") y <-
## Note: Mean no. genes in input = 1000
## Note: no. genes in output = 1000
## Note: estimated proportion of input genes in output = 1
head(y)
## rna k9a
## NR4A3 68.07237 10.7310010
## HSPA1B 47.19114 18.8135155
## DNAJB1 35.12799 2.4326983
## MIR133A1HG -27.36199 8.9061967
## HSPH1 25.83750 10.8922577
## CXCL2 24.76570 0.8459414
mitch can do unidimensional analysis if you provide it a single profile as a dataframe (not in a list).
mitch_import(rna,DEtype="edger") y <-
## The input is a single dataframe; one contrast only. Converting
## it to a list for you.
## Note: Mean no. genes in input = 1000
## Note: no. genes in output = 1000
## Note: estimated proportion of input genes in output = 1
head(y)
## x
## NR4A3 68.07237
## HSPA1B 47.19114
## DNAJB1 35.12799
## MIR133A1HG -27.36199
## HSPH1 25.83750
## CXCL2 24.76570
If the gene identifiers are not given in the rownames, then the column can be specified with the geneIDcol
parameter like this:
# first rearrange cols
rna
rna_mod <-$MyGeneIDs <- rownames(rna_mod)
rna_modrownames(rna_mod) <- seq(nrow(rna_mod))
head(rna_mod)
## logFC logCPM PValue adj.p.value MyGeneIDs
## 1 4.12734 5.09552 8.46507e-69 6.24341e-65 NR4A3
## 2 3.64685 7.42834 6.43968e-48 3.16639e-44 HSPA1B
## 3 2.35432 7.13208 7.44748e-36 2.74645e-32 DNAJB1
## 4 -1.02085 7.29935 4.34519e-28 1.28192e-24 MIR133A1HG
## 5 1.11729 6.03741 1.45377e-26 3.06351e-23 HSPH1
## 6 5.48158 2.88719 1.71515e-25 3.16252e-22 CXCL2
# now import with geneIDcol
mitch_import(rna_mod,DEtype="edgeR",geneIDcol="MyGeneIDs") y <-
## The input is a single dataframe; one contrast only. Converting
## it to a list for you.
## Note: Mean no. genes in input = 1000
## Note: no. genes in output = 1000
## Note: estimated proportion of input genes in output = 1
head(y)
## x
## NR4A3 68.07237
## HSPA1B 47.19114
## DNAJB1 35.12799
## MIR133A1HG -27.36199
## HSPH1 25.83750
## CXCL2 24.76570
By default, differential gene activity is scored using a supplied test statistic or directional p-value (D):
D = sgn(logFC) * -log10(p-value)
If this is not desired, then users can perform their own custom scoring procedure and import with DEtype="prescored"
.
There are many cases where the gene IDs don’t match the gene sets. To overcome this, mitch_import
also accepts a two-column table (gt here) that relates gene identifiers in the profiling data to those in the gene sets. In this example we can create some fake gene accession numbers to demonstrate this feature.
library("stringi")
# obtain vector of gene names
rownames(rna)
genenames <-# create fake accession numbers
paste("Gene0",stri_rand_strings(nrow(rna)*2, 6, pattern = "[0-9]"),sep="")
accessions <- head(unique(accessions),nrow(rna))
accessions <-# create a gene table file that relates gene names to accession numbers
data.frame(genenames,accessions)
gt <-
# now swap gene names for accessions
merge(rna,gt,by.x=0,by.y="genenames")
rna2 <-rownames(rna2) <- rna2$accessions
rna2[,2:5]
rna2 <-
merge(k9a,gt,by.x=0,by.y="genenames")
k9a2 <-rownames(k9a2) <- k9a2$accessions
k9a2[,2:5]
k9a2 <-
# now have a peek at the input data before importing
head(rna2,3)
## logFC logCPM PValue adj.p.value
## Gene0379337 0.296028 6.82814 3.84512e-04 1.46941e-02
## Gene0798939 -0.375440 4.71470 1.09120e-03 3.05432e-02
## Gene0145640 0.882624 8.12078 2.11945e-11 8.01642e-09
head(k9a2,3)
## logFC logCPM PValue adj.p.value
## Gene0379337 0.339535 3.67309 1.62925e-03 1.43363e-02
## Gene0798939 0.585837 3.66069 3.23724e-04 4.06552e-03
## Gene0145640 1.138700 2.78713 4.94270e-14 1.69263e-11
head(gt,3)
## genenames accessions
## 1 NR4A3 Gene0711975
## 2 HSPA1B Gene0926290
## 3 DNAJB1 Gene0542354
list("rna2"=rna2,"k9a2"=k9a2)
x <- mitch_import(x,DEtype="edgeR",geneTable=gt) y <-
## Note: Mean no. genes in input = 1000
## Note: no. genes in output = 1000
## Note: estimated proportion of input genes in output = 1
head(y,3)
## rna2 k9a2
## A2M 3.415090 2.788012
## AAAS -2.962096 3.489825
## ABRA 10.673777 13.306036
?mitch_import
provides more instructions on using this feature.
The mitch_calc
function performs multivariate enrichment analysis of the supplied gene sets in the scored profiling data. At its simpest form mitch_calc
function accepts the scored data as the first argument and the genesets as the second argument. Users can prioritise enrichments based on small adjusted p-values, by the observed effect size (“s.dist”, the magnitude of “s”, the enrichment score) or the standard deviation of the s scores. For most analyses, you will see more specific results when prioritising by effect, but it is worth looking at both.
Note that the number of parallel threads (“cores”) is set here to cores=2 but the default is to use just one. If you have many available threads, then a good number is 8-12. Note that the more CPU threads used, the more RAM is required, so monitor the available memory on your system.
# prioritisation by significance
mitch_calc(y,genesetsExample,priority="significance",cores=2) res <-
## Note: When prioritising by significance (ie: small
## p-values), large effect sizes might be missed.
# peek at the results
head(res$enrichment_result)
## set setSize
## 5 Biological oxidations 10
## 2 Antigen processing: Ubiquitination & Proteasome degradation 20
## 1 Adaptive Immune System 44
## 3 Asparagine N-linked glycosylation 17
## 4 Axon guidance 48
## pMANOVA s.rna2 s.k9a2 p.rna2 p.k9a2 s.dist
## 5 0.00899742 -0.555757576 -0.04141414 0.002429897 0.82165390 0.5572985
## 2 0.05970462 -0.055510204 0.28377551 0.670730417 0.02956405 0.2891538
## 1 0.52997174 0.000000000 0.09775580 1.000000000 0.27259775 0.0977558
## 3 0.53769811 -0.136137873 0.04589791 0.335578518 0.74549830 0.1436668
## 4 0.94522077 0.006959034 -0.02551646 0.935141485 0.76540436 0.0264484
## SD p.adjustMANOVA
## 5 0.36369573 0.0449871
## 2 0.23991123 0.1492616
## 1 0.06912379 0.6721226
## 3 0.12871874 0.6721226
## 4 0.02296364 0.9452208
# prioritisation by effect size
mitch_calc(y,genesetsExample,priority="effect",cores=2) res <-
## Note: Enrichments with large effect sizes may not be
## statistically significant.
head(res$enrichment_result)
## set setSize
## 5 Biological oxidations 10
## 2 Antigen processing: Ubiquitination & Proteasome degradation 20
## 3 Asparagine N-linked glycosylation 17
## 1 Adaptive Immune System 44
## 4 Axon guidance 48
## pMANOVA s.rna2 s.k9a2 p.rna2 p.k9a2 s.dist
## 5 0.00899742 -0.555757576 -0.04141414 0.002429897 0.82165390 0.5572985
## 2 0.05970462 -0.055510204 0.28377551 0.670730417 0.02956405 0.2891538
## 3 0.53769811 -0.136137873 0.04589791 0.335578518 0.74549830 0.1436668
## 1 0.52997174 0.000000000 0.09775580 1.000000000 0.27259775 0.0977558
## 4 0.94522077 0.006959034 -0.02551646 0.935141485 0.76540436 0.0264484
## SD p.adjustMANOVA
## 5 0.36369573 0.0449871
## 2 0.23991123 0.1492616
## 3 0.12871874 0.6721226
## 1 0.06912379 0.6721226
## 4 0.02296364 0.9452208
By default, gene sets with fewer than 10 members present in the profiling data are discarded. This threshold can be modified using the minsetsize option. There is no upper limit of gene set size.
mitch_calc(y,genesetsExample,priority="significance",minsetsize=5,cores=2) res <-
## Note: When prioritising by significance (ie: small
## p-values), large effect sizes might be missed.
By default, in downstream visualisation steps, charts are made from the top 50 gene sets, but this can be modified using the resrows option.
mitch_calc(y,genesetsExample,priority="significance",resrows=3,cores=2) res<-
## Note: When prioritising by significance (ie: small
## p-values), large effect sizes might be missed.
The HTML reports contain several plots as raster images and interactive charts which are useful as a first-pass visualisation. These can be generated like this:
mitch_report(res,"myreport.html",overwrite=TRUE)
## Dataset saved as " /tmp/RtmpufKX2y/myreport.rds ".
## processing file: mitch.Rmd
## output file: /tmp/RtmpLOtWVS/Rbuild2419b24b33b02/mitch/vignettes/mitch.knit.md
##
## Output created: /tmp/RtmpufKX2y/mitch_report.html
In case you want the charts in PDF format, for example for publications, these can be generated as such:
mitch_plots(res,outfile="mycharts.pdf")
For one-dimensional analysis, you can generate a network plot to visualise the similarity of the driver gene members for the top (20) up- and downregulated sets meeting a significance cut-off (FDR<0.05). These cut-offs can be altrered.
By default, this includes genes scored in the top and bottom tertiles, otherwise they are discarded. The intensity of the colour is proportional to the s.dist value (enrichment score). Circle size is proportional to the number of genes in the set. Line thickness is proportional to the Jaccard similarity value. This function works best after prioritisation with “effect” when running mitch_calc()
. Note that the circle (node) size and the line width shown in the legend is approximate, although the values shown are exactly the smallest and largest respectively. Note that this chart works best when the width is double the height, otherwise many of the long gene set names could be cut off. There is an element of stochasticity with regard to the network projection, so it could be a good idea to repeat it a few times until you get a nice layout.
Lists of enriched gene set members are shown.
networkplot(res,FDR=1,n_sets=10)
## Can't plot upregulated sets. Fewer than 5 found.
## Can't plot downregulated sets. Fewer than 5 found.
network_genes(res,FDR=1,n_sets=10)
## No significant upregulated sets to show.
## No significant upregulated sets to show.
## [[1]]
## NULL
sessionInfo()
## R Under development (unstable) (2025-03-13 r87965)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.2 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.21-bioc/R/lib/libRblas.so
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0 LAPACK version 3.12.0
##
## 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
##
## time zone: America/New_York
## tzcode source: system (glibc)
##
## attached base packages:
## [1] parallel stats4 stats graphics grDevices utils datasets
## [8] methods base
##
## other attached packages:
## [1] stringi_1.8.7
## [2] gtools_3.9.5
## [3] IlluminaHumanMethylationEPICanno.ilm10b4.hg19_0.6.0
## [4] IlluminaHumanMethylation450kanno.ilmn12.hg19_0.6.1
## [5] minfi_1.53.1
## [6] bumphunter_1.49.0
## [7] locfit_1.5-9.12
## [8] iterators_1.0.14
## [9] foreach_1.5.2
## [10] Biostrings_2.75.4
## [11] XVector_0.47.2
## [12] SummarizedExperiment_1.37.0
## [13] Biobase_2.67.0
## [14] MatrixGenerics_1.19.1
## [15] matrixStats_1.5.0
## [16] GenomicRanges_1.59.1
## [17] GenomeInfoDb_1.43.4
## [18] IRanges_2.41.3
## [19] S4Vectors_0.45.4
## [20] BiocGenerics_0.53.6
## [21] generics_0.1.3
## [22] HGNChelper_0.8.15
## [23] dplyr_1.1.4
## [24] mitch_1.19.5
##
## loaded via a namespace (and not attached):
## [1] splines_4.6.0 later_1.4.1
## [3] BiocIO_1.17.1 bitops_1.0-9
## [5] tibble_3.2.1 preprocessCore_1.69.0
## [7] XML_3.99-0.18 lifecycle_1.0.4
## [9] lattice_0.22-6 MASS_7.3-65
## [11] base64_2.0.2 scrime_1.3.5
## [13] magrittr_2.0.3 limma_3.63.11
## [15] sass_0.4.9 rmarkdown_2.29
## [17] jquerylib_0.1.4 yaml_2.3.10
## [19] httpuv_1.6.15 doRNG_1.8.6.1
## [21] askpass_1.2.1 DBI_1.2.3
## [23] RColorBrewer_1.1-3 abind_1.4-8
## [25] quadprog_1.5-8 purrr_1.0.4
## [27] RCurl_1.98-1.17 GenomeInfoDbData_1.2.14
## [29] rentrez_1.2.3 genefilter_1.89.0
## [31] annotate_1.85.0 svglite_2.1.3
## [33] DelayedMatrixStats_1.29.1 codetools_0.2-20
## [35] DelayedArray_0.33.6 xml2_1.3.8
## [37] tidyselect_1.2.1 farver_2.1.2
## [39] UCSC.utils_1.3.1 beanplot_1.3.1
## [41] illuminaio_0.49.0 GenomicAlignments_1.43.0
## [43] jsonlite_2.0.0 multtest_2.63.0
## [45] survival_3.8-3 systemfonts_1.2.1
## [47] tools_4.6.0 Rcpp_1.0.14
## [49] glue_1.8.0 gridExtra_2.3
## [51] SparseArray_1.7.7 xfun_0.51
## [53] HDF5Array_1.35.16 withr_3.0.2
## [55] fastmap_1.2.0 GGally_2.2.1
## [57] rhdf5filters_1.19.2 openssl_2.3.2
## [59] caTools_1.18.3 digest_0.6.37
## [61] R6_2.6.1 mime_0.13
## [63] colorspace_2.1-1 RSQLite_2.3.9
## [65] h5mread_0.99.4 tidyr_1.3.1
## [67] data.table_1.17.0 rtracklayer_1.67.1
## [69] httr_1.4.7 htmlwidgets_1.6.4
## [71] S4Arrays_1.7.3 ggstats_0.9.0
## [73] pkgconfig_2.0.3 gtable_0.3.6
## [75] blob_1.2.4 siggenes_1.81.0
## [77] htmltools_0.5.8.1 echarts4r_0.4.5
## [79] scales_1.3.0 kableExtra_1.4.0
## [81] png_0.1-8 knitr_1.50
## [83] rstudioapi_0.17.1 tzdb_0.5.0
## [85] reshape2_1.4.4 rjson_0.2.23
## [87] coda_0.19-4.1 statnet.common_4.11.0
## [89] nlme_3.1-168 curl_6.2.2
## [91] cachem_1.1.0 rhdf5_2.51.2
## [93] stringr_1.5.1 KernSmooth_2.23-26
## [95] AnnotationDbi_1.69.0 restfulr_0.0.15
## [97] GEOquery_2.75.0 pillar_1.10.1
## [99] grid_4.6.0 reshape_0.8.9
## [101] vctrs_0.6.5 gplots_3.2.0
## [103] promises_1.3.2 xtable_1.8-4
## [105] beeswarm_0.4.0 evaluate_1.0.3
## [107] readr_2.1.5 GenomicFeatures_1.59.1
## [109] cli_3.6.4 compiler_4.6.0
## [111] Rsamtools_2.23.1 rlang_1.1.5
## [113] crayon_1.5.3 rngtools_1.5.2
## [115] labeling_0.4.3 nor1mix_1.3-3
## [117] mclust_6.1.1 plyr_1.8.9
## [119] viridisLite_0.4.2 network_1.19.0
## [121] BiocParallel_1.41.5 munsell_0.5.1
## [123] Matrix_1.7-3 hms_1.1.3
## [125] sparseMatrixStats_1.19.0 bit64_4.6.0-1
## [127] ggplot2_3.5.1 Rhdf5lib_1.29.2
## [129] KEGGREST_1.47.0 statmod_1.5.0
## [131] shiny_1.10.0 memoise_2.0.1
## [133] bslib_0.9.0 bit_4.6.0
## [135] splitstackshape_1.4.8