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 Under development (unstable) (2024-10-21 r87258)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.1 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
#>
#> 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] stats graphics grDevices utils datasets methods base
#>
#> other attached packages:
#> [1] mitoClone2_1.13.0
#>
#> loaded via a namespace (and not attached):
#> [1] tidyselect_1.2.1 farver_2.1.2
#> [3] dplyr_1.1.4 blob_1.2.4
#> [5] Biostrings_2.75.0 bitops_1.0-9
#> [7] fastmap_1.2.0 RCurl_1.98-1.16
#> [9] VariantAnnotation_1.53.0 GenomicAlignments_1.43.0
#> [11] XML_3.99-0.17 digest_0.6.37
#> [13] lifecycle_1.0.4 KEGGREST_1.47.0
#> [15] RSQLite_2.3.7 magrittr_2.0.3
#> [17] compiler_4.5.0 rlang_1.1.4
#> [19] sass_0.4.9 tools_4.5.0
#> [21] utf8_1.2.4 yaml_2.3.10
#> [23] rtracklayer_1.67.0 Rhtslib_3.3.0
#> [25] knitr_1.48 S4Arrays_1.7.0
#> [27] bit_4.5.0 curl_5.2.3
#> [29] DelayedArray_0.33.0 RColorBrewer_1.1-3
#> [31] abind_1.4-8 BiocParallel_1.41.0
#> [33] BiocGenerics_0.53.0 grid_4.5.0
#> [35] stats4_4.5.0 fansi_1.0.6
#> [37] colorspace_2.1-1 ggplot2_3.5.1
#> [39] scales_1.3.0 SummarizedExperiment_1.37.0
#> [41] cli_3.6.3 rmarkdown_2.28
#> [43] crayon_1.5.3 generics_0.1.3
#> [45] httr_1.4.7 rjson_0.2.23
#> [47] DBI_1.2.3 cachem_1.1.0
#> [49] zlibbioc_1.53.0 splines_4.5.0
#> [51] parallel_4.5.0 AnnotationDbi_1.69.0
#> [53] formatR_1.14 XVector_0.47.0
#> [55] restfulr_0.0.15 matrixStats_1.4.1
#> [57] vctrs_0.6.5 Matrix_1.7-1
#> [59] jsonlite_1.8.9 VGAM_1.1-12
#> [61] IRanges_2.41.0 S4Vectors_0.45.0
#> [63] bit64_4.5.2 deepSNV_1.53.0
#> [65] GenomicFeatures_1.59.0 jquerylib_0.1.4
#> [67] glue_1.8.0 codetools_0.2-20
#> [69] gtable_0.3.6 GenomeInfoDb_1.43.0
#> [71] BiocIO_1.17.0 GenomicRanges_1.59.0
#> [73] UCSC.utils_1.3.0 munsell_0.5.1
#> [75] tibble_3.2.1 pillar_1.9.0
#> [77] htmltools_0.5.8.1 GenomeInfoDbData_1.2.13
#> [79] BSgenome_1.75.0 R6_2.5.1
#> [81] evaluate_1.0.1 lattice_0.22-6
#> [83] Biobase_2.67.0 highr_0.11
#> [85] png_0.1-8 Rsamtools_2.23.0
#> [87] pheatmap_1.0.12 memoise_2.0.1
#> [89] bslib_0.8.0 SparseArray_1.7.0
#> [91] xfun_0.48 MatrixGenerics_1.19.0
#> [93] pkgconfig_2.0.3