# IMPORTANT: this vignette can not be created if HiTC is not installed
if (!require("HiTC", quietly = TRUE)) {
  knitr::opts_chunk$set(eval = FALSE)
}
## 
## Attaching package: 'BiocGenerics'
## The following objects are masked from 'package:stats':
## 
##     IQR, mad, sd, var, xtabs
## The following objects are masked from 'package:base':
## 
##     anyDuplicated, append, as.data.frame, basename, cbind, colnames,
##     dirname, do.call, duplicated, eval, evalq, Filter, Find, get, grep,
##     grepl, intersect, is.unsorted, lapply, Map, mapply, match, mget,
##     order, paste, pmax, pmax.int, pmin, pmin.int, Position, rank,
##     rbind, Reduce, rownames, sapply, setdiff, sort, table, tapply,
##     union, unique, unsplit, which.max, which.min
## 
## Attaching package: 'S4Vectors'
## The following object is masked from 'package:utils':
## 
##     findMatches
## The following objects are masked from 'package:base':
## 
##     expand.grid, I, unname

Introduction

Hi-C is a sequencing-based molecular assay designed to measure intra and inter-chromosomal interactions between the DNA molecule. In particular, the identification of Topologically-Associated Domains (TADs), that is, of regions of the genome in which physical interactions are frequent, provides insight into the three-dimensional organization of a genome [2].

Hi-C data are in the form of two-dimensional contact maps, i.e., matrices whose \(i,j\) entry quantifies the intensity of the physical interaction between two genome regions \(i\) and \(j\) at the DNA level. In this vignette, we demonstrate the use of adjclust::hicClust to perform adjacency-constrained hierarchical agglomerative clustering (HAC) of Hi-C contact maps. The output of this function is a dendrogram, which can be cut to identify TADs. The algorithm used for adjacency-constrained (HAC) is described in [3,4].

library("adjclust")

Loading and displaying a sample Hi-C contact map

The data set hic_imr90_40_XX is an object of class HTCexp which has been obtained from the HiTC package [4]. It is a contact map corresponding to the first 500 x 500 bins on chromosome X vs chromosome X.

load(system.file("extdata", "hic_imr90_40_XX.rda", package = "adjclust"))

The script used to create this map can be found by executing the following command:

system.file("system/create_hic_chrXchrX.R", package="adjclust")
## [1] "/tmp/RtmpmMAagZ/Rinst101390451abf87/adjclust/system/create_hic_chrXchrX.R"

Now we have a look at the data.

HiTC::mapC(hic_imr90_40_XX)

Using hicClust

hicClust operates directly on objects of class HTCexp

fit <- hicClust(hic_imr90_40_XX)

It is also possible to work on binned data. Below we choose a bin size of \(5 \times 10^5\):

binned <- HiTC::binningC(hic_imr90_40_XX, binsize = 1e5)
fitB <- hicClust(binned)
fitB
## 
## Call:
## hicClust(binned)
## 
## Cluster method   : hicClust 
## Number of objects: 205
HiTC::mapC(binned)

The output is of class chac. In particular, it can be plotted as a dendrogram silently relying on the function plot.dendrogram:

plot(fitB, mode = "corrected")

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.

head(cbind(fitB$merge, fitB$gains))
##      [,1] [,2]
## [1,]   -3   -4
## [2,]   -2    1
## [3,]    2   -5
## [4,]   -1    3
## [5,]    4   -6
## [6,]  -17  -18

Other types of input

Contacts maps can also be stored as objects of class Matrix::dsCMatrix, or as plain text files. These types of input are also accepted as first argument to hicClust.

References

[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] Dixon J.R., et al (2012). Topological domains in mammalian genomes identified by analysis of chromatin interactions. Nature, 485(7398), 376.

[3] 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.

[4] Servant N., et al (2012). HiTC: Exploration of High-Throughput ‘C’ experiments. Bioinformatics, 28(21), 2843-2844.

Session information

sessionInfo()
## R version 4.3.2 (2023-10-31)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 22.04.3 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.20.so;  LAPACK version 3.10.0
## 
## locale:
##  [1] LC_CTYPE=fr_FR.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=fr_FR.UTF-8        LC_COLLATE=fr_FR.UTF-8    
##  [5] LC_MONETARY=fr_FR.UTF-8    LC_MESSAGES=fr_FR.UTF-8   
##  [7] LC_PAPER=fr_FR.UTF-8       LC_NAME=C                 
##  [9] LC_ADDRESS=C               LC_TELEPHONE=C            
## [11] LC_MEASUREMENT=fr_FR.UTF-8 LC_IDENTIFICATION=C       
## 
## time zone: Europe/Andorra
## tzcode source: system (glibc)
## 
## attached base packages:
## [1] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
## [1] adjclust_0.6.9       HiTC_1.40.0          GenomicRanges_1.48.0
## [4] GenomeInfoDb_1.32.4  IRanges_2.30.1       S4Vectors_0.34.0    
## [7] BiocGenerics_0.42.0 
## 
## loaded via a namespace (and not attached):
##  [1] SummarizedExperiment_1.26.1 gtable_0.3.4               
##  [3] capushe_1.1.2               rjson_0.2.21               
##  [5] xfun_0.40                   bslib_0.5.1                
##  [7] ggplot2_3.4.3               Biobase_2.56.0             
##  [9] lattice_0.20-45             vctrs_0.6.3                
## [11] tools_4.3.2                 bitops_1.0-7               
## [13] generics_0.1.3              parallel_4.3.2             
## [15] tibble_3.2.1                fansi_1.0.4                
## [17] highr_0.10                  pkgconfig_2.0.3            
## [19] Matrix_1.6-4                RColorBrewer_1.1-3         
## [21] sparseMatrixStats_1.8.0     lifecycle_1.0.3            
## [23] GenomeInfoDbData_1.2.8      compiler_4.3.2             
## [25] Rsamtools_2.12.0            Biostrings_2.64.1          
## [27] munsell_0.5.0               codetools_0.2-18           
## [29] htmltools_0.5.6             sass_0.4.7                 
## [31] RCurl_1.98-1.13             yaml_2.3.7                 
## [33] pillar_1.9.0                crayon_1.5.2               
## [35] jquerylib_0.1.4             MASS_7.3-57                
## [37] BiocParallel_1.30.4         DelayedArray_0.22.0        
## [39] cachem_1.0.8                viridis_0.6.4              
## [41] tidyselect_1.2.0            digest_0.6.33              
## [43] dplyr_1.1.2                 restfulr_0.0.15            
## [45] fastmap_1.1.1               grid_4.3.2                 
## [47] colorspace_2.1-0            cli_3.6.1                  
## [49] magrittr_2.0.3              XML_3.99-0.16              
## [51] utf8_1.2.3                  scales_1.2.1               
## [53] rmarkdown_2.24              XVector_0.36.0             
## [55] matrixStats_1.2.0           gridExtra_2.3              
## [57] evaluate_0.21               knitr_1.43                 
## [59] BiocIO_1.6.0                viridisLite_0.4.2          
## [61] rtracklayer_1.56.1          rlang_1.1.1                
## [63] dendextend_1.17.1           Rcpp_1.0.12                
## [65] glue_1.6.2                  rstudioapi_0.15.0          
## [67] jsonlite_1.8.7              R6_2.5.1                   
## [69] MatrixGenerics_1.8.1        GenomicAlignments_1.32.1   
## [71] zlibbioc_1.42.0