You should have identified true somatic variants in the mitochondrial (and nuclear) genome. The remaining vignettes of this package document how to get there. Here, we start with count matrices of the alternative and the reference alleles, across a number of sites of interest. Such data is available from two patients (P1, P2) as part of this package. Only P1 is analyzed as part of this vignette but the commented code for P2 should be fully functional.
load(system.file("data/M_P1.RData", package = "mitoClone2"))
load(system.file("data/N_P1.RData", package = "mitoClone2"))
P1 <- mutationCallsFromMatrix(as.matrix(M_P1), as.matrix(N_P1))
A first important step is to decide which mutations to include in the clustering. The default is to use all mutations that are covered in at least 20% of the cells, but this assignment can be changed manually.
The next step is to run SCITE or PhISCS to compute the most likely phylogenetic tree. PhISCS is bundled in this package, but the package needs to be run in an environment where gurobi and the gurobipy
python package are available. For example, you could set up a conda
environment that contains this package. SCITE does not require any additional licenses but may need to be manually compiled.
## this next step takes approx 4.1 minutes to run
tmpd <- tempdir()
dir.create(paste0(tmpd, "/p1"))
P1 <- varCluster(P1, tempfolder = paste0(tmpd, "/p1"), method = "SCITE")
#> Warning in max(which(x == 1)): no non-missing arguments to max; returning -Inf
This step can take a while to run. It computes a likely phylogenetic tree of all the mutations. In the case of SCITE, multiple equally likely tree can be produced. In it’s current state, this package simply selects the first tree from the list.
In many cases, the order of the leaves on these trees is arbitrary, because mutations systematically co-occur. We therefore cluster the mutations into clones. In detail, we take every every branch on the tree and then shuffle the order of mutations in that branch while re-calculating the likelihood. If swapping nodes leads to small changes in the likelihood, these nodes are then merged into a “clone”. The parameter min.lik
that controls the merging is set arbitrarily, see below for more information.
P1 <- clusterMetaclones(P1, min.lik = 1)
#> Warning in mutcalls@mut2clone[branches[[i]]] <-
#> as.integer(max(mutcalls@mut2clone) + : number of items to replace is not a
#> multiple of replacement length
This step also assigns each cell to the most likely clone, and provides an estimate of the likelihood. Refer to help(mutationCalls)
for more info on how these results are stored.
Finally, the clustering can be plotted.
The parameter min.lik
that controls the merging is set arbitrarily. In practice, the goal of these analyses is to group mutations into clones for subsequent analyses (such as differential expression analyses) and it may make sense to overwrite the result of clusterMetaclones
manually; for example, if a subclone defned on a mitochondrial mutation only should be treated as part of a more clearly defined upstream clone for differential expression analysis.
To overwrite the result of clusterMetaclones
, first retrieve the assignment of mutations to clones:
m2c <- mitoClone2:::getMut2Clone(P1)
print(m2c)
#> X1282GA X2537GA X3335TC X3350TC X5492TC X7527GDEL X8167TC X11196GA
#> 2 3 5 5 1 4 5 2
#> X11559GA X14462GA X16233AG EAPP SRSF2 CEBPA KLF7 root
#> 4 3 5 4 4 3 2 5
## To e.g. treat the mt:2537G>A and mt:14462:G>A mutations as a
## subclone distinct from CEBPA, we can assign them a new clonal
## identity
m2c[c("X2537GA", "X14462GA")] <- as.integer(6)
P1.new <- mitoClone2:::overwriteMetaclones(P1, m2c)
plotClones(P1.new)
#> R version 4.2.1 (2022-06-23)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 20.04.5 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.16-bioc/R/lib/libRblas.so
#> LAPACK: /home/biocbuild/bbs-3.16-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] stats graphics grDevices utils datasets methods base
#>
#> other attached packages:
#> [1] mitoClone2_1.4.0
#>
#> loaded via a namespace (and not attached):
#> [1] bitops_1.0-7 matrixStats_0.62.0
#> [3] bit64_4.0.5 RColorBrewer_1.1-3
#> [5] filelock_1.0.2 progress_1.2.2
#> [7] httr_1.4.4 GenomeInfoDb_1.34.0
#> [9] tools_4.2.1 bslib_0.4.0
#> [11] utf8_1.2.2 R6_2.5.1
#> [13] DBI_1.1.3 BiocGenerics_0.44.0
#> [15] colorspace_2.0-3 tidyselect_1.2.0
#> [17] prettyunits_1.1.1 bit_4.0.4
#> [19] curl_4.3.3 compiler_4.2.1
#> [21] cli_3.4.1 Biobase_2.58.0
#> [23] formatR_1.12 xml2_1.3.3
#> [25] DelayedArray_0.24.0 rtracklayer_1.58.0
#> [27] sass_0.4.2 scales_1.2.1
#> [29] deepSNV_1.44.0 rappdirs_0.3.3
#> [31] stringr_1.4.1 digest_0.6.30
#> [33] Rsamtools_2.14.0 rmarkdown_2.17
#> [35] XVector_0.38.0 pkgconfig_2.0.3
#> [37] htmltools_0.5.3 MatrixGenerics_1.10.0
#> [39] highr_0.9 dbplyr_2.2.1
#> [41] fastmap_1.1.0 BSgenome_1.66.0
#> [43] rlang_1.0.6 RSQLite_2.2.18
#> [45] VGAM_1.1-7 farver_2.1.1
#> [47] jquerylib_0.1.4 BiocIO_1.8.0
#> [49] generics_0.1.3 jsonlite_1.8.3
#> [51] BiocParallel_1.32.0 dplyr_1.0.10
#> [53] VariantAnnotation_1.44.0 RCurl_1.98-1.9
#> [55] magrittr_2.0.3 GenomeInfoDbData_1.2.9
#> [57] Matrix_1.5-1 Rcpp_1.0.9
#> [59] munsell_0.5.0 S4Vectors_0.36.0
#> [61] fansi_1.0.3 lifecycle_1.0.3
#> [63] stringi_1.7.8 yaml_2.3.6
#> [65] SummarizedExperiment_1.28.0 zlibbioc_1.44.0
#> [67] BiocFileCache_2.6.0 grid_4.2.1
#> [69] blob_1.2.3 parallel_4.2.1
#> [71] crayon_1.5.2 lattice_0.20-45
#> [73] Biostrings_2.66.0 splines_4.2.1
#> [75] GenomicFeatures_1.50.0 hms_1.1.2
#> [77] KEGGREST_1.38.0 knitr_1.40
#> [79] pillar_1.8.1 GenomicRanges_1.50.0
#> [81] rjson_0.2.21 codetools_0.2-18
#> [83] biomaRt_2.54.0 stats4_4.2.1
#> [85] XML_3.99-0.12 glue_1.6.2
#> [87] evaluate_0.17 png_0.1-7
#> [89] vctrs_0.5.0 gtable_0.3.1
#> [91] assertthat_0.2.1 cachem_1.0.6
#> [93] ggplot2_3.3.6 xfun_0.34
#> [95] restfulr_0.0.15 pheatmap_1.0.12
#> [97] tibble_3.1.8 GenomicAlignments_1.34.0
#> [99] AnnotationDbi_1.60.0 memoise_2.0.1
#> [101] IRanges_2.32.0 Rhtslib_2.0.0
#> [103] ellipsis_0.3.2