Compartmap: Direct inference of higher-order chromatin in single cells from scRNA-seq and scATAC-seq

Compartmap extends methods proposed by Fortin and Hansen 2015, Genome Biology (https://genomebiology.biomedcentral.com/articles/10.1186/s13059-015-0741-y) to perform direct inference of higher-order chromatin in single cells from scRNA-seq and scATAC-seq. Originally, Fortin and Hansen demonstrated that chromatin conformation could be inferred from (sc)ATAC-seq, bisulfite sequencing, DNase-seq and methylation arrays, similar to the results provided by HiC at the group level. Thus, in addition to the base information provided by the aforementioned assays, chromatin state could also be inferred.

Here, we propose a method to infer both group and single-cell level higher-order chromatin states from scRNA-seq and scATAC-seq. To accomplish this, we employ a James-Stein estimator (JSE) towards a global or targeted mean, using either chromsome or genome-wide information from scRNA-seq and scATAC-seq. Additionally, due to the sparsity of single-cell data, we employ a bootstrap procedure to quantify the uncertainty associated with the state and boundaries of inferred compartments. The output from compartmap can then be visualized directly, compared with orthogonal assay types, and/or embedded with something like UMAP or t-SNE. Further, to explore the higher-order interacting domains inferred from compartmap, we use a Random Matrix Theory (RMT) approach to resolve the “plaid-like” patterning, similar to what is observed in Hi-C and scHi-C.

Quick start with example data

Input

The expected input for compartmap is a RangedSummarizedExperiment object. These can be built using the built-in function importBigWig() if starting from BigWigs (recommended for scRNA-seq) or from a feature level object like a SingleCellExperiment with the rowRanges slot populated with the GRanges for each feature (see below in the examples).

Inferring higher-order chromatin domains at the group and single-cell level


#### Group level inference ####
#Process chr14 of the example K562 scRNA-seq data and infer higher-order chromatin at 1Mb resolution
k562_compartments <- scCompartments(k562_scrna_chr14,
                                    chr = "chr14",
                                    res = 1e6,
                                    group = TRUE,
                                    bootstrap = FALSE,
                                    genome = "hg19",
                                    assay = "rna")

#### Single-cell level inference ####
# To infer higher-order domains in single cells and quantifying sign coherence with the bootstrapping procedure, you can run:
# Sub-sample to 10 cells as an example
k562_scrna_chr14.sub <- k562_scrna_chr14[,sample(colnames(k562_scrna_chr14),
                                                 size = 10, replace = FALSE)]
k562_compartments.boot <- scCompartments(k562_scrna_chr14.sub,
                                         chr = "chr14",
                                         res = 1e6,
                                         group = FALSE,
                                         bootstrap = TRUE,
                                         num.bootstraps = 10,
                                         genome = "hg19",
                                         assay = "rna")

# Flip the domain sign if the sign coherence is discordant in 80% of the bootstraps
k562_compartments.boot.fix <- fixCompartments(k562_compartments.boot,
                                              min.conf = 0.8)

# Look at the first cell in the GRangesList object
k562_compartments.boot.fix[[1]]
## GRanges object with 89 ranges and 13 metadata columns:
##                             seqnames              ranges strand |        pc
##                                <Rle>           <IRanges>  <Rle> | <numeric>
##     chr14:19000000-19999999    chr14   19000000-19999999      * |  0.432049
##     chr14:20000000-20999999    chr14   20000000-20999999      * | -0.102012
##     chr14:21000000-21999999    chr14   21000000-21999999      * | -0.140159
##     chr14:22000000-22999999    chr14   22000000-22999999      * | -0.311821
##     chr14:23000000-23999999    chr14   23000000-23999999      * | -0.254601
##                         ...      ...                 ...    ... .       ...
##   chr14:103000000-103999999    chr14 103000000-103999999      * |  0.432049
##   chr14:104000000-104999999    chr14 104000000-104999999      * |  0.432049
##   chr14:105000000-105999999    chr14 105000000-105999999      * |  0.432049
##   chr14:106000000-106999999    chr14 106000000-106999999      * |  0.432049
##   chr14:107000000-107349539    chr14 107000000-107349539      * |  0.432049
##                             compartments     score boot.open boot.closed
##                              <character> <numeric> <numeric>   <numeric>
##     chr14:19000000-19999999         open  0.432049        10           0
##     chr14:20000000-20999999       closed -0.102012         4           6
##     chr14:21000000-21999999       closed -0.140159         4           6
##     chr14:22000000-22999999       closed -0.311821         0          10
##     chr14:23000000-23999999       closed -0.254601         0          10
##                         ...          ...       ...       ...         ...
##   chr14:103000000-103999999         open  0.432049        10           0
##   chr14:104000000-104999999         open  0.432049        10           0
##   chr14:105000000-105999999         open  0.432049        10           0
##   chr14:106000000-106999999         open  0.432049        10           0
##   chr14:107000000-107349539         open  0.432049        10           0
##                              conf.est conf.est.upperCI conf.est.lowerCI
##                             <numeric>        <numeric>        <numeric>
##     chr14:19000000-19999999  0.861234         1.000000         0.679113
##     chr14:20000000-20999999  0.572247         0.832889         0.311604
##     chr14:21000000-21999999  0.572247         0.832889         0.311604
##     chr14:22000000-22999999  0.861234         1.000000         0.679113
##     chr14:23000000-23999999  0.861234         1.000000         0.679113
##                         ...       ...              ...              ...
##   chr14:103000000-103999999  0.861234                1         0.679113
##   chr14:104000000-104999999  0.861234                1         0.679113
##   chr14:105000000-105999999  0.861234                1         0.679113
##   chr14:106000000-106999999  0.861234                1         0.679113
##   chr14:107000000-107349539  0.861234                1         0.679113
##                             flip.compartment flip.score flip.conf.est
##                                    <logical>  <numeric>     <numeric>
##     chr14:19000000-19999999            FALSE   0.432049      0.861234
##     chr14:20000000-20999999            FALSE  -0.102012      0.572247
##     chr14:21000000-21999999            FALSE  -0.140159      0.572247
##     chr14:22000000-22999999            FALSE  -0.311821      0.861234
##     chr14:23000000-23999999            FALSE  -0.254601      0.861234
##                         ...              ...        ...           ...
##   chr14:103000000-103999999            FALSE   0.432049      0.861234
##   chr14:104000000-104999999            FALSE   0.432049      0.861234
##   chr14:105000000-105999999            FALSE   0.432049      0.861234
##   chr14:106000000-106999999            FALSE   0.432049      0.861234
##   chr14:107000000-107349539            FALSE   0.432049      0.861234
##                             flip.conf.est.upperCI flip.conf.est.lowerCI
##                                         <numeric>             <numeric>
##     chr14:19000000-19999999              1.000000              0.679113
##     chr14:20000000-20999999              0.832889              0.311604
##     chr14:21000000-21999999              0.832889              0.311604
##     chr14:22000000-22999999              1.000000              0.679113
##     chr14:23000000-23999999              1.000000              0.679113
##                         ...                   ...                   ...
##   chr14:103000000-103999999                     1              0.679113
##   chr14:104000000-104999999                     1              0.679113
##   chr14:105000000-105999999                     1              0.679113
##   chr14:106000000-106999999                     1              0.679113
##   chr14:107000000-107349539                     1              0.679113
##   -------
##   seqinfo: 1 sequence from an unspecified genome; no seqlengths

Visualization of inferred chromatin domains

Once the data have been processed at either the group or single-cell level, one can visualize the results using the plotAB function in compartmap. Notably, we can include the confidence intervals and median, chromosome-wide confidence estimate derived from the bootstrap procedure for sign coherence. At 50%, this suggests that estimates are evenly split between open and closed states. This may be due to data sparsity or heterogeneity in the data. One possible approach to resolve this is to increase the number of bootstraps performed if initially set low (e.g. 10). Alternatively, it may be a region that is worth investigating for your data set.

Importing bigWigs as input to compartmap

The currently recommended input files to compartmap for scRNA-seq are single-cell bigWigs, though does work with a feature/counts based object as demonstrated in the next section. Single-cell bigWigs can be generated through several tools, such as deeptools (https://deeptools.readthedocs.io/en/latest/). To import bigWigs, we can use the importBigWig function in compartmap. This will read in a bigWig file and optionally summarize to an arbitrary bin size. The bin size used in the compartmap manuscript was 1kb and is what we do here as well.

Starting with a feature or counts-based object

In the cases where we do not have or can’t start with bigWigs (e.g. scATAC), we can start with a feature-level or counts object (e.g. SingleCellExperiment). The two things that must be there are making sure rowRanges and colnames are set for each feature and cell/sample. We will use the scATAC-seq from K562 as an example of how the object should look, but will also show one way to add rowRanges to a SingleCellExperiment, which works the same way for a SummarizedExperiment object.

Showing the rowRanges slot is a GRanges for each feature

But if we don’t have rowRanges for something like a SingleCellExperiment we are working with, we need to generate them. Thus, we will show an example of how to add rowRanges from a GTF file to a SingleCellExperiment.

Higher-order chromatin interaction maps

Another interesting aspect we can derive from scRNA and scATAC is the higher-order interacting domains through denoising of the correlation matrices using a Random Matrix Theory approach. This is often represented with the “plaid-like” patterning shown in Hi-C and scHi-C approaches where stronger correlations (e.g. greater intensity of red) indicates interacting domains relative to lesser correlation. We can do something similar here.