The goal of fastreeR is to provide functions for calculating distance matrix, building phylogenetic tree or performing hierarchical clustering between samples, directly from a VCF or FASTA file.
To install fastreeR
package:
You should allocate minimum 10kb per sample per variant of RAM for
the JVM. The more RAM you allocate, the faster the execution will be
(less pauses for garbage collection). In order to allocate RAM, a
special parameter needs to be passed while JVM initializes. JVM
parameters can be passed by setting java.parameters
option.
The -Xmx
parameter, followed (without space) by an integer
value and a letter, is used to tell JVM what is the maximum amount of
heap RAM that it can use. The letter in the parameter (uppercase or
lowercase), indicates RAM units. For example, parameters
-Xmx1024m
or -Xmx1024M
or -Xmx1g
or -Xmx1G
, all allocate 1 Gigabyte or 1024 Megabytes of
maximum RAM for JVM.
We download, in a temporary location, a small vcf file from 1K
project, with around 150 samples and 100k variants (SNPs and INDELs). We
use BiocFileCache
for this retrieval process so that it is
not repeated needlessly. If for any reason we cannot download, we use
the small sample vcf from fastreeR
package.
bfc <- BiocFileCache::BiocFileCache(ask = FALSE)
tempVcfUrl <-
paste0("https://ftp.1000genomes.ebi.ac.uk/vol1/ftp/data_collections/",
"1000_genomes_project/release/20190312_biallelic_SNV_and_INDEL/",
"supporting/related_samples/",
"ALL.chrX.shapeit2_integrated_snvindels_v2a_related_samples_27022019.",
"GRCh38.phased.vcf.gz")
tempVcf <- BiocFileCache::bfcquery(bfc,field = "rname", "tempVcf")$rpath[1]
if(is.na(tempVcf)) {
tryCatch(
{ tempVcf <- BiocFileCache::bfcadd(bfc,"tempVcf",fpath=tempVcfUrl)[[1]]
},
error=function(cond) {
tempVcf <- system.file("extdata", "samples.vcf.gz", package="fastreeR")
},
warning=function(cond) {
tempVcf <- system.file("extdata", "samples.vcf.gz", package="fastreeR")
}
)
}
if(file.size(tempVcf) == 0L) {
tempVcf <- system.file("extdata", "samples.vcf.gz", package="fastreeR")
}
We download, in temporary location, some small bacterial genomes. We
use BiocFileCache
for this retrieval process so that it is
not repeated needlessly. If for any reason we cannot download, we use
the small sample fasta from fastreeR
package.
tempFastasUrls <- c(
#Mycobacterium liflandii
paste0("https://ftp.ncbi.nih.gov/genomes/refseq/bacteria/",
"Mycobacterium_liflandii/latest_assembly_versions/",
"GCF_000026445.2_ASM2644v2/GCF_000026445.2_ASM2644v2_genomic.fna.gz"),
#Pelobacter propionicus
paste0("https://ftp.ncbi.nih.gov/genomes/refseq/bacteria/",
"Pelobacter_propionicus/latest_assembly_versions/",
"GCF_000015045.1_ASM1504v1/GCF_000015045.1_ASM1504v1_genomic.fna.gz"),
#Rickettsia prowazekii
paste0("https://ftp.ncbi.nih.gov/genomes/refseq/bacteria/",
"Rickettsia_prowazekii/latest_assembly_versions/",
"GCF_000022785.1_ASM2278v1/GCF_000022785.1_ASM2278v1_genomic.fna.gz"),
#Salmonella enterica
paste0("https://ftp.ncbi.nih.gov/genomes/refseq/bacteria/",
"Salmonella_enterica/reference/",
"GCF_000006945.2_ASM694v2/GCF_000006945.2_ASM694v2_genomic.fna.gz"),
#Staphylococcus aureus
paste0("https://ftp.ncbi.nih.gov/genomes/refseq/bacteria/",
"Staphylococcus_aureus/reference/",
"GCF_000013425.1_ASM1342v1/GCF_000013425.1_ASM1342v1_genomic.fna.gz")
)
tempFastas <- list()
for (i in seq(1,5)) {
tempFastas[[i]] <- BiocFileCache::bfcquery(bfc,field = "rname",
paste0("temp_fasta",i))$rpath[1]
if(is.na(tempFastas[[i]])) {
tryCatch(
{ tempFastas[[i]] <-
BiocFileCache::bfcadd(bfc, paste0("temp_fasta",i),
fpath=tempFastasUrls[i])[[1]]
},
error=function(cond) {
tempFastas <- system.file("extdata", "samples.fasta.gz",
package="fastreeR")
break
},
warning=function(cond) {
tempFastas <- system.file("extdata", "samples.fasta.gz",
package="fastreeR")
break
}
)
}
if(!file.exists(tempFastas[[i]])) {
tempFastas <- system.file("extdata", "samples.fasta.gz",
package="fastreeR")
break
}
if(file.size(tempFastas[[i]]) == 0L) {
tempFastas <- system.file("extdata", "samples.fasta.gz",
package="fastreeR")
break
}
}
The most time consuming process is calculating distances between samples. Assign more processors in order to speed up this operation.
We note two distinct groups of distances. One around of distance value 0.05 and the second around distance value 0.065.
fastreeR::dist2tree
Notice that the generated tree is ultrametric.
myVcfTree <- fastreeR::dist2tree(inputDist = myVcfDist)
plot(ape::read.tree(text = myVcfTree), direction = "down", cex = 0.3)
ape::add.scale.bar()
ape::axisPhylo(side = 2)
Of course the same can be achieved directly from the vcf file, without calculating distances.
myVcfTree <- fastreeR::vcf2tree(inputFile = tempVcf, threads = 2)
plot(ape::read.tree(text = myVcfTree), direction = "down", cex = 0.3)
ape::add.scale.bar()
ape::axisPhylo(side = 2)
As expected from the histogram of distances, two groups of samples also emerge in the tree. The two branches, one at height around 0.055 and the second around height 0.065, are clearly visible.
stats::hclust
For comparison, we generate a tree by using stats
package and distances calculated by fastreeR
.
Although it does not initially look very similar, because it is not ultrametric, it is indeed quite the same tree. We note again the two groups (two branches) of samples and the 4 samples, possibly clones, that they show very close distances between them.
We can identify the two groups of samples, apparent from the
hierarchical tree, by using dist2clusters
or
vcf2clusters
or tree2clusters
. By playing a
little with the cutHeight
parameter, we find that a value
of cutHeight=0.067
cuts the tree into two branches. The
first group contains 106 samples and the second 44.
myVcfClust <- fastreeR::dist2clusters(inputDist = myVcfDist, cutHeight = 0.067)
#> ..done.
if (length(myVcfClust) > 1) {
tree <- myVcfClust[[1]]
clusters <- myVcfClust[[2]]
tree
clusters
}
#> [1] "1 106 HG00124 HG00153 HG00247 HG00418 HG00427 HG00501 HG00512 HG00577 HG00578 HG00635 HG00702 HG00716 HG00733 HG00866 HG00983 HG01195 HG01274 HG01278 HG01322 HG01347 HG01452 HG01453 HG01473 HG01477 HG01480 HG01482 HG01483 HG01590 HG01983 HG01995 HG02024 HG02046 HG02218 HG02288 HG02344 HG02347 HG02363 HG02372 HG02377 HG02381 HG02387 HG02388 HG02524 HG02525 HG02781 HG03487 HG03606 HG03618 HG03621 HG03633 HG03639 HG03650 HG03656 HG03699 HG03700 HG03715 HG03723 HG03761 HG03794 HG03797 HG03799 HG03806 HG03811 HG03842 HG03845 HG03847 HG03901 HG03904 HG03929 HG03948 HG03972 HG03982 HG03988 HG04024 HG04037 HG04050 HG04053 HG04055 HG04058 HG04114 HG04127 HG04128 HG04132 HG04135 HG04147 HG04149 HG04150 HG04174 HG04191 HG04192 NA07346 NA11993 NA12891 NA12892 NA19660 NA19675 NA19685 NA19737 NA19797 NA19798 NA20336 NA20344 NA20526 NA20871 NA20893 NA20898"
#> [2] "2 44 HG02478 HG02762 HG02869 HG02964 HG02965 HG03033 HG03034 HG03076 HG03249 HG03250 HG03306 HG03307 HG03309 HG03312 HG03339 HG03361 HG03373 HG03383 HG03408 HG03454 HG03493 HG03508 HG03566 HG03569 HG03574 HG03582 NA18487 NA19150 NA19240 NA19311 NA19313 NA19373 NA19381 NA19382 NA19396 NA19444 NA19453 NA19469 NA19470 NA19985 NA20313 NA20322 NA20341 NA20361"
Similar analysis we can perform when we have samples represented as sequences in a fasta file.
Use of the downloaded sample fasta file :
Or use the provided by fastreeR
fasta file of 48
bacterial RefSeq :
utils::sessionInfo()
#> R version 4.4.2 (2024-10-31)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.1 LTS
#>
#> Matrix products: default
#> BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
#> LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
#>
#> 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=en_US.UTF-8
#> [9] LC_ADDRESS=en_US.UTF-8 LC_TELEPHONE=en_US.UTF-8
#> [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=en_US.UTF-8
#>
#> time zone: Etc/UTC
#> tzcode source: system (glibc)
#>
#> attached base packages:
#> [1] grid stats graphics grDevices utils datasets methods
#> [8] base
#>
#> other attached packages:
#> [1] BiocFileCache_2.15.0 dbplyr_2.5.0 ape_5.8
#> [4] fastreeR_1.11.0 BiocStyle_2.35.0
#>
#> loaded via a namespace (and not attached):
#> [1] sass_0.4.9 utf8_1.2.4 generics_0.1.3
#> [4] stringi_1.8.4 RSQLite_2.3.8 lattice_0.22-6
#> [7] digest_0.6.37 magrittr_2.0.3 evaluate_1.0.1
#> [10] dynamicTreeCut_1.63-1 fastmap_1.2.0 blob_1.2.4
#> [13] R.oo_1.27.0 jsonlite_1.8.9 R.utils_2.12.3
#> [16] DBI_1.2.3 BiocManager_1.30.25 httr_1.4.7
#> [19] purrr_1.0.2 fansi_1.0.6 jquerylib_0.1.4
#> [22] cli_3.6.3 rlang_1.1.4 R.methodsS3_1.8.2
#> [25] bit64_4.5.2 withr_3.0.2 cachem_1.1.0
#> [28] yaml_2.3.10 tools_4.4.2 parallel_4.4.2
#> [31] memoise_2.0.1 dplyr_1.1.4 filelock_1.0.3
#> [34] curl_6.0.1 buildtools_1.0.0 vctrs_0.6.5
#> [37] R6_2.5.1 lifecycle_1.0.4 stringr_1.5.1
#> [40] bit_4.5.0 pkgconfig_2.0.3 rJava_1.0-11
#> [43] pillar_1.9.0 bslib_0.8.0 glue_1.8.0
#> [46] Rcpp_1.0.13-1 xfun_0.49 tibble_3.2.1
#> [49] tidyselect_1.2.1 sys_3.4.3 knitr_1.49
#> [52] htmltools_0.5.8.1 nlme_3.1-166 rmarkdown_2.29
#> [55] maketools_1.3.1 compiler_4.4.2