The hardware and bandwidth for this mirror is donated by dogado GmbH, the Webhosting and Full Service-Cloud Provider. Check out our Wordpress Tutorial.
If you wish to report a bug, or if you are interested in having us mirror your free-software or open-source project, please feel free to contact us at mirror[@]dogado.de.

Introduction

hicream is an R package designed to perform Hi-c data differential analysis. More specifically, it performs a pixel-level differential analysis using diffHic and, using a two-dimensionnal connectivity constrained clustering, renders a post hoc bound on True Discovery Proportion for each cluster. This method allows to identify differential genomic regions and quantifies signal in those regions.

library("hicream")
## Loading required package: reticulate
# checking if Python and Python modules are available to avoid errors in vignette building
modules_avail <- reticulate::py_available(initialize = TRUE) &&
  reticulate::py_module_available("sklearn") &&
  reticulate::py_module_available("kneebow") &&
  reticulate::py_module_available("pandas") &&
  reticulate::py_module_available("numpy")

How to load external data?

Using hicream requires to load data that corresponds to Hi-C matrices and their index, the latter in bed format. In the following code, the paths to Hi-C matrices and the index file path are used in the loadData function alongside the chromosome number. The option normalize = TRUE allows to perform a cyclic LOESS normalization. The output is an InteractionSet object.

replicates <- 1:2
cond <- "90"
allBegins <- interaction(expand.grid(replicates, cond), sep = "-")
allBegins <- as.character(allBegins)
chromosome <- 1
nbChr <- 1
allMat <- sapply(allBegins, function(ab) {
  matFile <- paste0("Rep", ab, "-chr", chromosome, "_200000.bed")
})
index <- system.file("extdata", "index.200000.longest18chr.abs.bed",
                     package = "hicream")
format <- rep("HiC-Pro", length(replicates) * length(cond) * nbChr)
binsize <- 200000
files <- system.file("extdata", unlist(allMat), package = "hicream")
exData <- loadData(files, index, chromosome, normalize = TRUE)
##  chr(0)
## 
##  chr(0)
## 

pighic dataset

The pighic dataset has been produced using 6 Hi-C matrices (3 in each condition) obtained from two different developmental stages of pig embryos (90 and 110 days of gestation) corresponding to chromosome 1. The dataset contains two elements, the first is an output of the loadData function corresponding to normalized data. The second is the vector giving, for each matrix, the condition it belongs to.

data("pighic")
head(pighic)
## $data
## class: InteractionSet 
## dim: 21 6 
## metadata(1): width
## assays(2): counts offset
## rownames: NULL
## rowData names(2): anchor1.id anchor2.id
## colnames: NULL
## colData names(1): totals
## type: ReverseStrictGInteractions
## regions: 6
## 
## $conditions
## [1] "90"  "90"  "90"  "110" "110" "110"

Perform hicream test

Once data have been loaded, the three steps of the analysis can be performed.

Use performTest function

In this first step, the output of a loadData object is used, with a vector indicating the condition of each sample. Here, we use data stored in pighic. performTest outputs the result of the diffHic pixel-level differential analysis.

resdiff <- performTest(pighic$data, pighic$conditions)
resdiff
## Testing for differential interactions using diffHic.
## 
##   region1 region2   p.value     p.adj       logFC
## 1    1125    1125 0.7318482 0.8538229 -0.01830721
## 2    1125    1126 0.0562111 0.2342558 -0.16338780
##  [ reached 'max' / getOption("max.print") -- omitted 19 rows ]
summary(resdiff)
## Summary of diffHic results.
##     p.value              p.adj              logFC          
##  Min.   :0.0002064   Min.   :0.004334   Min.   :-0.163388  
##  1st Qu.:0.0669302   1st Qu.:0.234256   1st Qu.:-0.034500  
##  Median :0.2751093   Median :0.525209   Median : 0.003099  
##  Mean   :0.3960997   Mean   :0.580792   Mean   : 0.022271  
##  3rd Qu.:0.6924986   3rd Qu.:0.853823   3rd Qu.: 0.040991  
##  Max.   :0.9537715   Max.   :0.953771   Max.   : 0.349596

Several plotting options allow to look at the \(p\)-value map, adjusted \(p\)-value map or log-fold-change map.

plot(resdiff)

plot(resdiff, which_plot = "p.adj")

plot(resdiff, which_plot = "logFC")

Perform two-dimensionnal connectivity-constrained clustering

The AggloClust2D function uses the output of loadData function in order to build a connectivity graph of pixels and perform a two-dimensionnal hierarchical clustering under the connectivity constraint. The function renders a clustering corresponding to the optimal number of clusters found by the elbow heuristic. However, a clustering corresponding to any other number of clusters (chosen by the user) can be obtained by specifying a value for the input argument nbClust.

if (modules_avail) {
  res2D <- AggloClust2D(pighic$data)
  res2D 
  summary(res2D)
}
## Summary of 2D constrained clustering results.
##             Length Class  Mode     
## merge       40     -none- numeric  
## height      20     -none- numeric  
## order       21     -none- numeric  
## labels      21     -none- numeric  
## method       1     -none- character
## call         2     -none- call     
## dist.method  1     -none- character

Plotting the output shows the tree that corresponds to the hierarchy of clusters.

if (modules_avail) {
  plot(res2D)
}

Perform post hoc inference

Post hoc inference is performed by the postHoc function, using the output of performTest, and a clustering (the latter can be obtained with AggloClust2D). The user sets a level of test confidence \(\alpha\), typically equal to \(0.05\).

if (modules_avail) {
  clusters <- res2D$clustering
  alpha <- 0.05
  resposthoc <- postHoc(resdiff, clusters, alpha)
  resposthoc
  summary(resposthoc)
}
## Summary of post hoc results.
##      TPRate           p.value              p.adj              logFC          
##  Min.   :0.00000   Min.   :0.0002064   Min.   :0.004334   Min.   :-0.163388  
##  1st Qu.:0.00000   1st Qu.:0.0669302   1st Qu.:0.234256   1st Qu.:-0.034500  
##  Median :0.00000   Median :0.2751093   Median :0.525209   Median : 0.003099  
##  Mean   :0.04762   Mean   :0.3960997   Mean   :0.580792   Mean   : 0.022271  
##  3rd Qu.:0.16667   3rd Qu.:0.6924986   3rd Qu.:0.853823   3rd Qu.: 0.040991  
##  Max.   :0.16667   Max.   :0.9537715   Max.   :0.953771   Max.   : 0.349596  
##   propPoslogFC   
##  Min.   :0.0000  
##  1st Qu.:0.5000  
##  Median :0.5000  
##  Mean   :0.5238  
##  3rd Qu.:0.6667  
##  Max.   :1.0000

Plotting the output of postHoc function shows the minimal proportion of True Discoveries in each cluster.

if (modules_avail) plot(resposthoc)

Session Information

sessionInfo()
## R version 4.5.1 (2025-06-13)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.2 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=fr_FR.UTF-8        LC_COLLATE=C              
##  [5] LC_MONETARY=fr_FR.UTF-8    LC_MESSAGES=en_US.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/Paris
## tzcode source: system (glibc)
## 
## attached base packages:
## [1] stats     graphics  grDevices utils     datasets  methods   base     
## 
## other attached packages:
## [1] hicream_0.0.1     reticulate_1.42.0
## 
## loaded via a namespace (and not attached):
##   [1] tidyselect_1.2.1            viridisLite_0.4.2          
##   [3] dplyr_1.1.4                 farver_2.1.2               
##   [5] viridis_0.6.5               Biostrings_2.76.0          
##   [7] bitops_1.0-9                fastmap_1.2.0              
##   [9] RCurl_1.98-1.17             GenomicAlignments_1.44.0   
##  [11] XML_3.99-0.18               digest_0.6.37              
##  [13] lifecycle_1.0.4             capushe_1.1.2              
##  [15] statmod_1.5.0               magrittr_2.0.3             
##  [17] compiler_4.5.1              rlang_1.1.6                
##  [19] sass_0.4.10                 tools_4.5.1                
##  [21] yaml_2.3.10                 rtracklayer_1.68.0         
##  [23] knitr_1.50                  Rhtslib_3.4.0              
##  [25] labeling_0.4.3              S4Arrays_1.8.0             
##  [27] curl_6.4.0                  here_1.0.1                 
##  [29] DelayedArray_0.34.1         plyr_1.8.9                 
##  [31] RColorBrewer_1.1-3          abind_1.4-8                
##  [33] BiocParallel_1.42.1         withr_3.0.2                
##  [35] BiocGenerics_0.54.0         grid_4.5.1                 
##  [37] stats4_4.5.1                Rhdf5lib_1.30.0            
##  [39] edgeR_4.6.2                 ggplot2_3.5.2              
##  [41] scales_1.4.0                MASS_7.3-65                
##  [43] SummarizedExperiment_1.38.1 cli_3.6.5                  
##  [45] rmarkdown_2.29              crayon_1.5.3               
##  [47] auk_0.8.2                   generics_0.1.4             
##  [49] metapod_1.16.0              reshape2_1.4.4             
##  [51] rjson_0.2.23                httr_1.4.7                 
##  [53] rhdf5_2.52.1                cachem_1.1.0               
##  [55] stringr_1.5.1               splines_4.5.1              
##  [57] parallel_4.5.1              XVector_0.48.0             
##  [59] restfulr_0.0.16             matrixStats_1.5.0          
##  [61] vctrs_0.6.5                 Matrix_1.7-3               
##  [63] jsonlite_2.0.0              IRanges_2.42.0             
##  [65] S4Vectors_0.46.0            dendextend_1.19.0          
##  [67] adjclust_0.6.10             locfit_1.5-9.12            
##  [69] limma_3.64.1                jquerylib_0.1.4            
##  [71] glue_1.8.0                  codetools_0.2-20           
##  [73] stringi_1.8.7               gtable_0.3.6               
##  [75] GenomeInfoDb_1.44.0         GenomicRanges_1.60.0       
##  [77] BiocIO_1.18.0               UCSC.utils_1.4.0           
##  [79] tibble_3.3.0                pillar_1.10.2              
##  [81] rappdirs_0.3.3              rhdf5filters_1.20.0        
##  [83] csaw_1.42.0                 htmltools_0.5.8.1          
##  [85] GenomeInfoDbData_1.2.14     BSgenome_1.76.0            
##  [87] R6_2.6.1                    sparseMatrixStats_1.20.0   
##  [89] diffHic_1.40.0              rprojroot_2.0.4            
##  [91] evaluate_1.0.4              lattice_0.22-7             
##  [93] Biobase_2.68.0              png_0.1-8                  
##  [95] Rsamtools_2.24.0            bslib_0.9.0                
##  [97] Rcpp_1.0.14                 InteractionSet_1.36.1      
##  [99] gridExtra_2.3               SparseArray_1.8.0          
## [101] xfun_0.52                   MatrixGenerics_1.20.0      
## [103] pkgconfig_2.0.3

These binaries (installable software) and packages are in development.
They may not be fully stable and should be used with caution. We make no claims about them.
Health stats visible at Monitor.