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.
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:
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.
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
):
library(nullranges)
set.seed(123)
mgi <- matchRanges(focal = binPairs[binPairs$looped],
pool = binPairs[!binPairs$looped],
covar = ~ctcfSignal + distance + n_sites,
method = 'stratified')
mgi
## 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:
library(plyranges)
library(ggplot2)
## Summarize ctcfSignal by n_sites
mgi %>%
regions() %>%
group_by(n_sites) %>%
summarize(ctcfSignal = mean(ctcfSignal)) %>%
as.data.frame() %>%
ggplot(aes(x = n_sites, y = ctcfSignal)) +
geom_line() +
geom_point(shape = 21, stroke = 1, fill = 'white') +
theme_minimal() +
theme(panel.border = element_rect(color = 'black',
fill = NA))
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
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:
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:
## Generate a randomly selected set from all binPairs
all <- c(focal(mgi), pool(mgi))
set.seed(123)
random <- all[sample(1:length(all), length(mgi), replace = FALSE)]
## Calculate the percent of convergent CTCF sites for each group
g1 <- (sum(focal(mgi)$convergent) / length(focal(mgi))) * 100
g2 <- (sum(pool(mgi)$convergent) / length(pool(mgi))) * 100
g3 <- (sum(random$convergent) / length(random)) * 100
g4 <- (sum(mgi$convergent) / length(mgi)) * 100
## Visualize
barplot(height = c(g1, g2, g3, g4),
names.arg = c('looped\n(focal)', 'unlooped\n(pool)',
'all pairs\n(random)', 'unlooped\n(matched)'),
col = c('#1F78B4', '#33A02C', 'orange2', '#A6CEE3'),
ylab = "Convergent CTCF Sites (%)",
main = "Testing the Convergence Rule",
border = NA,
las = 1)
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.
## R version 4.1.1 (2021-08-10)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.3 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.14-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.14-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.1 plyranges_1.14.0
## [3] nullrangesData_1.0.0 ExperimentHub_2.2.0
## [5] AnnotationHub_3.2.0 BiocFileCache_2.2.0
## [7] dbplyr_2.1.1 ggplot2_3.3.5
## [9] plotgardener_1.0.1 nullranges_1.0.1
## [11] InteractionSet_1.22.0 SummarizedExperiment_1.24.0
## [13] Biobase_2.54.0 MatrixGenerics_1.6.0
## [15] matrixStats_0.61.0 GenomicRanges_1.46.0
## [17] GenomeInfoDb_1.30.0 IRanges_2.28.0
## [19] S4Vectors_0.32.2 BiocGenerics_0.40.0
##
## loaded via a namespace (and not attached):
## [1] colorspace_2.0-2 rjson_0.2.20
## [3] ellipsis_0.3.2 ggridges_0.5.3
## [5] mclust_5.4.8 XVector_0.34.0
## [7] farver_2.1.0 bit64_4.0.5
## [9] interactiveDisplayBase_1.32.0 AnnotationDbi_1.56.1
## [11] fansi_0.5.0 mvtnorm_1.1-3
## [13] cachem_1.0.6 knitr_1.36
## [15] jsonlite_1.7.2 speedglm_0.3-3
## [17] Rsamtools_2.10.0 png_0.1-7
## [19] shiny_1.7.1 BiocManager_1.30.16
## [21] compiler_4.1.1 httr_1.4.2
## [23] assertthat_0.2.1 Matrix_1.3-4
## [25] fastmap_1.1.0 later_1.3.0
## [27] prettyunits_1.1.1 htmltools_0.5.2
## [29] tools_4.1.1 gtable_0.3.0
## [31] glue_1.5.0 GenomeInfoDbData_1.2.7
## [33] dplyr_1.0.7 rappdirs_0.3.3
## [35] Rcpp_1.0.7 jquerylib_0.1.4
## [37] vctrs_0.3.8 Biostrings_2.62.0
## [39] strawr_0.0.9 rtracklayer_1.54.0
## [41] xfun_0.28 stringr_1.4.0
## [43] mime_0.12 lifecycle_1.0.1
## [45] restfulr_0.0.13 XML_3.99-0.8
## [47] zlibbioc_1.40.0 MASS_7.3-54
## [49] scales_1.1.1 hms_1.1.1
## [51] promises_1.2.0.1 parallel_4.1.1
## [53] RColorBrewer_1.1-2 yaml_2.2.1
## [55] curl_4.3.2 memoise_2.0.0
## [57] yulab.utils_0.0.4 sass_0.4.0
## [59] stringi_1.7.5 RSQLite_2.2.8
## [61] BiocVersion_3.14.0 highr_0.9
## [63] BiocIO_1.4.0 filelock_1.0.2
## [65] BiocParallel_1.28.0 rlang_0.4.12
## [67] pkgconfig_2.0.3 bitops_1.0-7
## [69] pracma_2.3.3 evaluate_0.14
## [71] lattice_0.20-45 purrr_0.3.4
## [73] GenomicAlignments_1.30.0 ks_1.13.2
## [75] labeling_0.4.2 bit_4.0.4
## [77] tidyselect_1.1.1 plyr_1.8.6
## [79] magrittr_2.0.1 R6_2.5.1
## [81] generics_0.1.1 DelayedArray_0.20.0
## [83] DBI_1.1.1 pillar_1.6.4
## [85] withr_2.4.2 KEGGREST_1.34.0
## [87] RCurl_1.98-1.5 tibble_3.1.6
## [89] crayon_1.4.2 KernSmooth_2.23-20
## [91] utf8_1.2.2 rmarkdown_2.11
## [93] progress_1.2.2 data.table_1.14.2
## [95] blob_1.2.2 digest_0.6.28
## [97] xtable_1.8-4 httpuv_1.6.3
## [99] gridGraphics_0.5-1 munsell_0.5.0
## [101] ggplotify_0.1.0 bslib_0.3.1