In this vignette, we demonstrate the generation of covariate-matched null ranges by using the matchRanges() function to test the “covergence rule” of CTCF-bound chromatin loops, first described in Rao et al. 2014.

Background

In the human genome, chromatin loops play an important role in regulating gene expression by connecting regulatory loci to gene promoters. The anchors of these loops are bound by the protein CTCF, which is required for loop formation and maintenance. CTCF binds a specific non-palindromic DNA sequence motif. And the direction of this motif at loop anchors is non-random. That is, CTCF binding motifs found at the anchors of loops tend to be oriented towards the center of the loop more often then would be expected by chance. This phenomenon, first described by Rao et al in 2014, is known as the “convergence rule”. It was initially discovered by first filtering for loops that had a single CTCF binding site at each end of the loop. Here we use matchRanges to reinvestigate the convergence using all loops.

We will use hg19_10kb_ctcfBoundBinPairs data from the nullrangesData package which contains features from the GM12878 cell line aligned to hg19 as well as CTCF motif data from the CTCF data package available on Bioconductor. hg19_10kb_ctcfBoundBinPairs is a GInteractions object with all pairwise combinations of CTCF-bound 10Kb bins within 1Mb annotated with the following features:

  • The total CTCF signal in each bin.
  • The number of CTCF sites in each bin.
  • The distance between bin pairs.
  • Whether at least one CTCF site is convergent between each bin pair.
  • The presence or absence of a loop between each bin pair.

Using these annotations and the matchRanges() function, we can compare CTCF motif orientations between pairs of genomic regions that are 1) connected by loops, 2) not connected by loops, 3) randomly chosen, or 4) not connected by loops, but matched for the strength of CTCF sites and distance between loop anchors.

Matching with matchRanges()

Before we generate our null ranges, let’s take a look at our example dataset:

## StrictGInteractions object with 198120 interactions and 5 metadata columns:
##            seqnames1             ranges1     seqnames2             ranges2 |
##                <Rle>           <IRanges>         <Rle>           <IRanges> |
##        [1]      chr1       230001-240000 ---      chr1       520001-530000 |
##        [2]      chr1       230001-240000 ---      chr1       710001-720000 |
##        [3]      chr1       230001-240000 ---      chr1       800001-810000 |
##        [4]      chr1       230001-240000 ---      chr1       840001-850000 |
##        [5]      chr1       230001-240000 ---      chr1       870001-880000 |
##        ...       ...                 ... ...       ...                 ... .
##   [198116]      chrX 154310001-154320000 ---      chrX 154370001-154380000 |
##   [198117]      chrX 154310001-154320000 ---      chrX 155250001-155260000 |
##   [198118]      chrX 154320001-154330000 ---      chrX 154370001-154380000 |
##   [198119]      chrX 154320001-154330000 ---      chrX 155250001-155260000 |
##   [198120]      chrX 154370001-154380000 ---      chrX 155250001-155260000 |
##               looped ctcfSignal  n_sites  distance convergent
##            <logical>  <numeric> <factor> <integer>  <logical>
##        [1]     FALSE    5.18038        2    290000      FALSE
##        [2]     FALSE    5.46775        2    480000       TRUE
##        [3]     FALSE    7.30942        2    570000      FALSE
##        [4]     FALSE    7.34338        2    610000      FALSE
##        [5]     FALSE    6.31338        3    640000       TRUE
##        ...       ...        ...      ...       ...        ...
##   [198116]     FALSE    6.79246        2     60000      FALSE
##   [198117]     FALSE    6.12447        3    940000       TRUE
##   [198118]     FALSE    7.40868        2     50000       TRUE
##   [198119]     FALSE    7.00936        3    930000      FALSE
##   [198120]     FALSE    6.73402        3    880000       TRUE
##   -------
##   regions: 20612 ranges and 5 metadata columns
##   seqinfo: 23 sequences from hg19 genome

Let’s start by defining our focal set (i.e. looped bin-pairs), our pool set (i.e un-looped bin-pairs), and our covariates of interest (i.e. ctcfSignal and distance):

## MatchedGInteractions object with 3104 interactions and 5 metadata columns:
##          seqnames1             ranges1     seqnames2             ranges2 |
##              <Rle>           <IRanges>         <Rle>           <IRanges> |
##      [1]     chr11   62160001-62170000 ---     chr11   62190001-62200000 |
##      [2]     chr17     7890001-7900000 ---     chr17     7990001-8000000 |
##      [3]     chr22   36460001-36470000 ---     chr22   36680001-36690000 |
##      [4]     chr11     1560001-1570000 ---     chr11     1710001-1720000 |
##      [5]     chr12   31780001-31790000 ---     chr12   31830001-31840000 |
##      ...       ...                 ... ...       ...                 ... .
##   [3100]      chr2   20980001-20990000 ---      chr2   21180001-21190000 |
##   [3101]     chr19   14310001-14320000 ---     chr19   14540001-14550000 |
##   [3102]      chr2   20060001-20070000 ---      chr2   20840001-20850000 |
##   [3103]     chr14 100750001-100760000 ---     chr14 100850001-100860000 |
##   [3104]     chr11     1840001-1850000 ---     chr11     2170001-2180000 |
##             looped ctcfSignal  n_sites  distance convergent
##          <logical>  <numeric> <factor> <integer>  <logical>
##      [1]     FALSE    8.12860        2     30000      FALSE
##      [2]     FALSE    8.72976        3    100000       TRUE
##      [3]     FALSE    9.51410        4    220000       TRUE
##      [4]     FALSE    8.64759        3    150000      FALSE
##      [5]     FALSE    8.43412        2     50000       TRUE
##      ...       ...        ...      ...       ...        ...
##   [3100]     FALSE    8.10867        2    200000       TRUE
##   [3101]     FALSE    8.13800        2    230000      FALSE
##   [3102]     FALSE    8.46727        2    780000      FALSE
##   [3103]     FALSE    8.52685        4    100000       TRUE
##   [3104]     FALSE    8.39022        3    330000       TRUE
##   -------
##   regions: 20612 ranges and 5 metadata columns
##   seqinfo: 23 sequences from hg19 genome

When the focal and pool arguments are GInteractions objects, matchRanges() returns a MatchedGInteractions object. The MatchedGInteractions class extends GInteractions so all of the same operations can be applied:

Assessing quality of matching

We can get a quick summary of the matching quality with overview():

## MatchedGInteractions object: 
##        set      N ctcfSignal.mean ctcfSignal.sd distance.mean distance.sd
##      focal   3104             8.3          0.67        320000      230000
##    matched   3104             8.3          0.71        320000      250000
##       pool 195016             7.9          0.85        490000      290000
##  unmatched 191912             7.8          0.85        490000      290000
##  n_sites.0 n_sites.1 n_sites.2 n_sites.3 n_sites.4 n_sites.5 n_sites.6
##          0         0      2167       734       164        32         4
##          0         0      2191       719       157        30         3
##          0         0    152318     34992      5971       944       158
##          0         0    150127     34273      5814       914       155
##  n_sites.7+ ps.mean ps.sd
##           3   0.026 0.016
##           4   0.026 0.016
##         633   0.016 0.013
##         629   0.015 0.013
## --------
## focal - matched: 
##  ctcfSignal.mean ctcfSignal.sd distance.mean distance.sd n_sites.0 n_sites.1
##         -2.8e-05        -0.047          -250      -26000         0         0
##  n_sites.2 n_sites.3 n_sites.4 n_sites.5 n_sites.6 n_sites.7+  ps.mean    ps.sd
##        -24        15         7         2         1         -1 -1.7e-07 -7.4e-07

In addition to providing a printed overview, the overview data can be extracted for convenience. For example, the quality property shows the absolute value of the mean difference between focal and matched sets. Therefore, the lower this value, the better the matching quality:

## [1] 1.7e-07

Visualizing matching results

Let’s visualize overall matching quality by plotting propensity scores for the focal, pool, and matched sets:

Log transformations can be applied to ‘x’, ‘y’, or both (c('x', 'y')) for plotting functions to make it easier to assess quality. It is clear that the matched set is very well matched to the focal set:

We can ensure that covariate distributions have been matched appropriately by using the covariates() function to extract matched covariates along with patchwork and plotCovarite to visualize all distributions:

Compare CTCF site orientation

Using our matched ranges, we can compare the percent of looped pairs with at least one convergent CTCF site against unlooped pairs, randomly selected pairs, and pairs that are unlooped but have been matched for our covariates. The accessor function focal() and pool() can be used to conveniently extract these matched sets:

As shown above the convergent rule holds true even when controlling for CTCF signal strength and bin pair distance. Greater than 90% of looped pairs contain convergent CTCF sites, compared to only ~55% among non-looped subsets.

Session information

## R version 4.2.1 (2022-06-23)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.5 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.16-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.16-bioc/R/lib/libRlapack.so
## 
## 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       
## 
## attached base packages:
## [1] grid      stats4    stats     graphics  grDevices utils     datasets 
## [8] methods   base     
## 
## other attached packages:
##  [1] patchwork_1.1.2             plyranges_1.18.0           
##  [3] nullrangesData_1.3.0        ExperimentHub_2.6.0        
##  [5] AnnotationHub_3.6.0         BiocFileCache_2.6.0        
##  [7] dbplyr_2.2.1                ggplot2_3.3.6              
##  [9] plotgardener_1.4.0          nullranges_1.4.0           
## [11] InteractionSet_1.26.0       SummarizedExperiment_1.28.0
## [13] Biobase_2.58.0              MatrixGenerics_1.10.0      
## [15] matrixStats_0.62.0          GenomicRanges_1.50.0       
## [17] GenomeInfoDb_1.34.0         IRanges_2.32.0             
## [19] S4Vectors_0.36.0            BiocGenerics_0.44.0        
## 
## loaded via a namespace (and not attached):
##   [1] colorspace_2.0-3              rjson_0.2.21                 
##   [3] ellipsis_0.3.2                ggridges_0.5.4               
##   [5] mclust_6.0.0                  XVector_0.38.0               
##   [7] farver_2.1.1                  bit64_4.0.5                  
##   [9] interactiveDisplayBase_1.36.0 AnnotationDbi_1.60.0         
##  [11] fansi_1.0.3                   mvtnorm_1.1-3                
##  [13] codetools_0.2-18              cachem_1.0.6                 
##  [15] knitr_1.40                    jsonlite_1.8.3               
##  [17] speedglm_0.3-4                Rsamtools_2.14.0             
##  [19] png_0.1-7                     shiny_1.7.3                  
##  [21] BiocManager_1.30.19           compiler_4.2.1               
##  [23] httr_1.4.4                    assertthat_0.2.1             
##  [25] Matrix_1.5-1                  fastmap_1.1.0                
##  [27] cli_3.4.1                     later_1.3.0                  
##  [29] prettyunits_1.1.1             htmltools_0.5.3              
##  [31] tools_4.2.1                   gtable_0.3.1                 
##  [33] glue_1.6.2                    GenomeInfoDbData_1.2.9       
##  [35] dplyr_1.0.10                  rappdirs_0.3.3               
##  [37] Rcpp_1.0.9                    jquerylib_0.1.4              
##  [39] vctrs_0.5.0                   Biostrings_2.66.0            
##  [41] strawr_0.0.9                  rtracklayer_1.58.0           
##  [43] xfun_0.34                     stringr_1.4.1                
##  [45] mime_0.12                     lifecycle_1.0.3              
##  [47] restfulr_0.0.15               XML_3.99-0.12                
##  [49] zlibbioc_1.44.0               MASS_7.3-58.1                
##  [51] scales_1.2.1                  hms_1.1.2                    
##  [53] promises_1.2.0.1              parallel_4.2.1               
##  [55] RColorBrewer_1.1-3            yaml_2.3.6                   
##  [57] curl_4.3.3                    memoise_2.0.1                
##  [59] yulab.utils_0.0.5             sass_0.4.2                   
##  [61] stringi_1.7.8                 RSQLite_2.2.18               
##  [63] BiocVersion_3.16.0            highr_0.9                    
##  [65] BiocIO_1.8.0                  filelock_1.0.2               
##  [67] BiocParallel_1.32.0           rlang_1.0.6                  
##  [69] pkgconfig_2.0.3               bitops_1.0-7                 
##  [71] pracma_2.4.2                  evaluate_0.17                
##  [73] lattice_0.20-45               purrr_0.3.5                  
##  [75] GenomicAlignments_1.34.0      ks_1.13.5                    
##  [77] labeling_0.4.2                bit_4.0.4                    
##  [79] tidyselect_1.2.0              magrittr_2.0.3               
##  [81] R6_2.5.1                      generics_0.1.3               
##  [83] DelayedArray_0.24.0           DBI_1.1.3                    
##  [85] pillar_1.8.1                  withr_2.5.0                  
##  [87] KEGGREST_1.38.0               RCurl_1.98-1.9               
##  [89] tibble_3.1.8                  crayon_1.5.2                 
##  [91] KernSmooth_2.23-20            utf8_1.2.2                   
##  [93] rmarkdown_2.17                progress_1.2.2               
##  [95] data.table_1.14.4             blob_1.2.3                   
##  [97] digest_0.6.30                 xtable_1.8-4                 
##  [99] httpuv_1.6.6                  gridGraphics_0.5-1           
## [101] munsell_0.5.0                 ggplotify_0.1.0              
## [103] bslib_0.4.0