# IMPORTANT: this vignette is not created if snpStats is not installed
if (!require("snpStats")) {
knitr::opts_chunk$set(eval = FALSE)
}
## Loading required package: snpStats
## Loading required package: survival
## Loading required package: Matrix
##
## Attaching package: 'Matrix'
## The following object is masked from 'package:S4Vectors':
##
## expand
In this vignette we demonstrate the use of snpClust
function in the adjclust
package. snpClust
performs adjacency-constrained hierarchical clustering of single
nucleotide polymorphisms (SNPs), where the similarity between SNPs is
defined by linkage disequilibrium (LD).
This function implements the algorithm described in [1]. It is an extension of the algorithm described in [3,4]. Denoting by p the number of SNPs to cluster and assuming that the similarity between SNPs whose indices are more distant than h, its time complexity is O(p(log (p) + h)), and its space complexity is O(hp).
The beginning of this vignette closely follows the “LD vignette” of the SnpStats package [2]. First, we load genotype data.
We focus on the ceph.1mb
data.
geno <- ceph.1mb[, -316] ## drop one SNP leading to one missing LD value
p <- ncol(geno)
nSamples <- nrow(geno)
geno
## A SnpMatrix with 90 rows and 602 columns
## Row names: NA06985 ... NA12892
## Col names: rs5993821 ... rs5747302
These data are drawn from the International HapMap Project and concern 602 SNPs1 over a 1Mb region of chromosome 22 in sample of 90 Europeans.
We can compute and display the LD between these SNPs.
The snpClust
function can handle genotype data as an
input:
## Warning in run.snpClust(x, h = h, stats = stats): Forcing the LD similarity to
## be smaller than or equal to 1
## Note: 135 merges with non increasing heights.
Note that due to numerical errors in the LD estimation, some of the estimated LD values may be slightly larger than 1. These values are rounded to 1 internally.
The above figure suggests that the LD signal is concentrated close to
the diagonal. We can focus on a diagonal band with the bandwidth
parameter h
:
## Warning in run.snpClust(x, h = h, stats = stats): Forcing the LD similarity to
## be smaller than or equal to 1
## Note: 133 merges with non increasing heights.
##
## Call:
## snpClust(geno, h = 100, stats = "R.squared")
##
## Cluster method : snpClust
## Number of objects: 602
The output of the snpClust
is of class
chac
. In particular, it can be plotted as a dendrogram
silently relying on the function plot.dendrogram
:
## Warning:
## Detected reversals in dendrogram: mode = 'corrected', 'within-disp' or 'total-disp' might be more relevant.
Moreover, the output contains an element named merge
which describes the successive merges of the clustering, and an element
gains
which gives the improvement in the criterion
optimized by the clustering at each successive merge.
## [,1] [,2]
## [1,] -1 -2
## [2,] -255 -256
## [3,] -488 -489
## [4,] -487 3
## [5,] -486 4
## [6,] -234 -235
In this section we show how the snpClust
function can
also be applied directly to LD values.
h <- 100
ld.ceph <- snpStats::ld(geno, stats = "R.squared", depth = h, symmetric = TRUE)
image(ld.ceph, lwd = 0)
Note that we have forced the snpStats::ld
function to
return a symmetric matrix. We can apply snpClust
directly
to this LD matrix (of class Matrix::dsCMatrix
):
## Note: forcing the diagonal of the LD similarity matrix to be 1
## Warning in run.snpClust(x, h = h, stats = stats): Forcing the LD similarity to
## be smaller than or equal to 1
## Note: 133 merges with non increasing heights.
snpClust
also handles inputs of class
base::matrix
:
## Note: 133 merges with non increasing heights.
[1] Ambroise C., Dehman A., Neuvial P., Rigaill G., and Vialaneix N. (2019). Adjacency-constrained hierarchical clustering of a band similarity matrix with application to genomics. Algorithms for Molecular Biology, 14, 22.
[2] Clayton D. (2015). snpStats: SnpMatrix and XSnpMatrix classes and methods. R package version 1.20.0
[3] Dehman A., Ambroise C., Neuvial P. (2015). Performance of a blockwise approach in variable selection using linkage disequilibrium information. BMC Bioinformatics, 16, 148.
[4] Randriamihamison N., Vialaneix N., and Neuvial P. (2021). Applicability and interpretability of Ward’s hierarchical agglomerative clustering with or without contiguity constraints. Journal of Classification, 38, 363–389.
## 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=C
## [9] LC_ADDRESS=C LC_TELEPHONE=C
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## time zone: Etc/UTC
## tzcode source: system (glibc)
##
## attached base packages:
## [1] stats4 stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] snpStats_1.57.0 Matrix_1.7-1 survival_3.7-0
## [4] adjclust_0.6.10 HiTC_1.51.0 GenomicRanges_1.59.0
## [7] GenomeInfoDb_1.43.0 IRanges_2.41.0 S4Vectors_0.45.0
## [10] BiocGenerics_0.53.1 generics_0.1.3
##
## loaded via a namespace (and not attached):
## [1] SummarizedExperiment_1.37.0 capushe_1.1.2
## [3] gtable_0.3.6 rjson_0.2.23
## [5] xfun_0.49 bslib_0.8.0
## [7] ggplot2_3.5.1 Biobase_2.67.0
## [9] lattice_0.22-6 vctrs_0.6.5
## [11] tools_4.4.2 bitops_1.0-9
## [13] curl_6.0.0 parallel_4.4.2
## [15] fansi_1.0.6 tibble_3.2.1
## [17] highr_0.11 pkgconfig_2.0.3
## [19] RColorBrewer_1.1-3 sparseMatrixStats_1.19.0
## [21] lifecycle_1.0.4 GenomeInfoDbData_1.2.13
## [23] compiler_4.4.2 Rsamtools_2.23.0
## [25] Biostrings_2.75.1 munsell_0.5.1
## [27] codetools_0.2-20 htmltools_0.5.8.1
## [29] sys_3.4.3 buildtools_1.0.0
## [31] sass_0.4.9 RCurl_1.98-1.16
## [33] yaml_2.3.10 pillar_1.9.0
## [35] crayon_1.5.3 jquerylib_0.1.4
## [37] MASS_7.3-61 BiocParallel_1.41.0
## [39] cachem_1.1.0 DelayedArray_0.33.1
## [41] viridis_0.6.5 abind_1.4-8
## [43] digest_0.6.37 restfulr_0.0.15
## [45] splines_4.4.2 maketools_1.3.1
## [47] fastmap_1.2.0 grid_4.4.2
## [49] colorspace_2.1-1 cli_3.6.3
## [51] SparseArray_1.7.1 magrittr_2.0.3
## [53] S4Arrays_1.7.1 utf8_1.2.4
## [55] XML_3.99-0.17 UCSC.utils_1.3.0
## [57] scales_1.3.0 rmarkdown_2.29
## [59] XVector_0.47.0 httr_1.4.7
## [61] matrixStats_1.4.1 gridExtra_2.3
## [63] evaluate_1.0.1 knitr_1.48
## [65] BiocIO_1.17.0 viridisLite_0.4.2
## [67] rtracklayer_1.67.0 rlang_1.1.4
## [69] dendextend_1.18.1 Rcpp_1.0.13-1
## [71] glue_1.8.0 jsonlite_1.8.9
## [73] R6_2.5.1 MatrixGenerics_1.19.0
## [75] GenomicAlignments_1.43.0 zlibbioc_1.52.0
We have dropped SNP rs2401075 because it produced a missing value due to the lack of genetic diversity in the considered sample.↩︎