Contents

1 Introduction

Gene regulatory networks model the underlying gene regulation hierarchies that drive gene expression and cell states. The main functions of this package are to construct gene regulatory networks and infer transcription factor (TF) activity at the single cell level by integrating scATAC-seq and scRNA-seq data and incorporating of public bulk TF ChIP-seq data.

There are three related packages: epiregulon, epiregulon.extra and epiregulon.archr, the two of which are available through Bioconductor and the last of which is only available through github. The basic epiregulon package takes in gene expression and chromatin accessibility matrices in the form of SingleCellExperiment objects, constructs gene regulatory networks (“regulons”) and outputs the activity of transcription factors at the single cell level. The epiregulon.extra package provides a suite of tools for enrichment analysis of target genes, visualization of target genes and transcription factor activity, and network analysis which can be run on the epiregulon output. If the users would like to start from ArchR projects instead of SingleCellExperiment objects, they may choose to use epiregulon.archr package, which allows for seamless integration with the ArchR package, and continue with the rest of the workflow offered in epiregulon.extra.

For full documentation of the epiregulon package, please refer to the epiregulon book.

2 Installation

if (!require("BiocManager", quietly = TRUE))
    install.packages("BiocManager")
 
BiocManager::install("epiregulon")

Alternatively, you could install from github

devtools::install_github(repo ='xiaosaiyao/epiregulon')

Load package

library(epiregulon)

3 Data preparation

Prior to using epiregulon, single cell preprocessing needs to performed by user’s favorite methods. The following components are required:
1. Peak matrix from scATAC-seq containing the chromatin accessibility information
2. Gene expression matrix from either paired or unpaired scRNA-seq. RNA-seq integration needs to be performed for unpaired dataset.
3. Dimensionality reduction matrix from either single modality dataset or joint scRNA-seq and scATAC-seq

This tutorial demonstrates the basic functions of epiregulon, using the reprogram-seq dataset which can be downloaded from the scMultiome package. In this example, prostate cancer cells (LNCaP) were infected in separate wells with viruses encoding 4 transcription factors (NKX2-1, GATA6, FOXA1 and FOXA2) and a positive control (mNeonGreen) before pooling. The identity of the infected transcription factors was tracked through cell hashing (available in the field hash_assignment of the colData) and serves as the ground truth.

# load the MAE object
library(scMultiome)
mae <- scMultiome::reprogramSeq()

# extract peak matrix
PeakMatrix <- mae[["PeakMatrix"]]

# extract expression matrix
GeneExpressionMatrix <- mae[["GeneExpressionMatrix"]]
rownames(GeneExpressionMatrix) <- rowData(GeneExpressionMatrix)$name

# define the order of hash_assigment
GeneExpressionMatrix$hash_assignment <- 
  factor(as.character(GeneExpressionMatrix$hash_assignment),
         levels = c("HTO10_GATA6_UTR", "HTO2_GATA6_v2", 
                    "HTO8_NKX2.1_UTR", "HTO3_NKX2.1_v2",
                    "HTO1_FOXA2_v2", "HTO4_mFOXA1_v2", "HTO6_hFOXA1_UTR", 
                    "HTO5_NeonG_v2"))
# extract dimensional reduction matrix
reducedDimMatrix <- reducedDim(mae[['TileMatrix500']], "LSI_ATAC")

# transfer UMAP_combined from TileMatrix to GeneExpressionMatrix
reducedDim(GeneExpressionMatrix, "UMAP_Combined") <- 
  reducedDim(mae[['TileMatrix500']], "UMAP_Combined")

Visualize singleCellExperiment by UMAP


scater::plotReducedDim(GeneExpressionMatrix, 
                       dimred = "UMAP_Combined", 
                       text_by = "Clusters", 
                       colour_by = "Clusters")

4 Quick start

4.1 Retrieve bulk TF ChIP-seq binding sites

First, we retrieve a GRangesList object containing the binding sites of all the transcription factors and co-regulators. These binding sites are derived from bulk ChIP-seq data in the ChIP-Atlas and ENCODE databases. For the same transcription factor, multiple ChIP-seq files from different cell lines or tissues are merged. For further information on how these peaks are derived, please refer to ?epiregulon::getTFMotifInfo. Currently, human genomes hg19 and hg38 and mouse mm10 are supported.

grl <- getTFMotifInfo(genome = "hg38")
#> Retrieving chip-seq data, version 2
#> see ?scMultiome and browseVignettes('scMultiome') for documentation
#> loading from cache
grl
#> GRangesList object of length 1377:
#> $AEBP2
#> GRanges object with 2700 ranges and 0 metadata columns:
#>          seqnames            ranges strand
#>             <Rle>         <IRanges>  <Rle>
#>      [1]     chr1        9792-10446      *
#>      [2]     chr1     942105-942400      *
#>      [3]     chr1     984486-984781      *
#>      [4]     chr1   3068932-3069282      *
#>      [5]     chr1   3069411-3069950      *
#>      ...      ...               ...    ...
#>   [2696]     chrY   8465261-8465730      *
#>   [2697]     chrY 11721744-11722260      *
#>   [2698]     chrY 11747448-11747964      *
#>   [2699]     chrY 19302661-19303134      *
#>   [2700]     chrY 19985662-19985982      *
#>   -------
#>   seqinfo: 25 sequences from an unspecified genome; no seqlengths
#> 
#> ...
#> <1376 more elements>

We recommend the use of ChIP-seq data over motif for estimating TF occupancy. However, if the user would like to start from motifs, it is possible to switch to the motif mode. The user can provide the peak regions as a GRanges object and the location of the motifs will be annotated based on Cisbp from chromVARmotifs (human_pwms_v2 and mouse_pwms_v2, version 0.2)

grl.motif <- getTFMotifInfo(genome = "hg38",  
                            mode = "motif", 
                            peaks = rowRanges(PeakMatrix))

4.3 Add TF motif binding to peaks

The next step is to add the TF binding information by overlapping regions of the peak matrix with the bulk chip-seq database. The output is a data frame object with three columns:

  1. idxATAC - index of the peak in the peak matrix
  2. idxTF - index in the gene expression matrix corresponding to the transcription factor
  3. tf - name of the transcription factor

overlap <- addTFMotifInfo(grl = grl, 
                          p2g = p2g, 
                          peakMatrix = PeakMatrix)
#> Computing overlap...
#> Success!
head(overlap)
#>   idxATAC idxTF    tf
#> 1       1    16  ATF1
#> 2       1    17  ATF2
#> 3       1    18  ATF3
#> 4       1    21  ATF7
#> 5       1    35 BRCA2
#> 6       1    36  BRD4

4.4 Generate regulons

A DataFrame, representing the inferred regulons, is then generated. The DataFrame consists of ten columns:

  1. idxATAC - index of the peak in the peak matrix
  2. chr - chromosome number
  3. start - start position of the peak
  4. end - end position of the peak
  5. idxRNA - index in the gene expression matrix corresponding to the target gene
  6. target - name of the target gene
  7. distance - distance between the transcription start site of the target gene and the middle of the peak
  8. idxTF - index in the gene expression matrix corresponding to the transcription factor
  9. tf - name of the transcription factor
  10. corr - correlation between target gene expression and the chromatin accessibility at the peak. if cluster labels are provided, this field is a matrix with columns names corresponding to correlation across all cells and for each of the clusters.

regulon <- getRegulon(p2g = p2g, 
                      overlap = overlap)
regulon
#> DataFrame with 5214284 rows and 12 columns
#>           idxATAC         chr     start       end    idxRNA  target      corr
#>         <integer> <character> <integer> <integer> <integer> <array>  <matrix>
#> 1               1        chr1    817121    817621        24  SAMD11 -0.388417
#> 2               1        chr1    817121    817621        24  SAMD11 -0.388417
#> 3               1        chr1    817121    817621        24  SAMD11 -0.388417
#> 4               1        chr1    817121    817621        24  SAMD11 -0.388417
#> 5               1        chr1    817121    817621        24  SAMD11 -0.388417
#> ...           ...         ...       ...       ...       ...     ...       ...
#> 5214280    126586        chrX 155071227 155071727     36420  FUNDC2  0.560009
#> 5214281    126586        chrX 155071227 155071727     36420  FUNDC2  0.560009
#> 5214282    126590        chrX 155228844 155229344     36426   CLIC2  0.787924
#> 5214283    126590        chrX 155228844 155229344     36426   CLIC2  0.787924
#> 5214284    126590        chrX 155228844 155229344     36426   CLIC2  0.787924
#>         p_val_peak_gene FDR_peak_gene  distance     idxTF          tf
#>                <matrix>      <matrix> <integer> <integer> <character>
#> 1             0.0391887        0.9055    108108        16        ATF1
#> 2             0.0391887        0.9055    108108        17        ATF2
#> 3             0.0391887        0.9055    108108        18        ATF3
#> 4             0.0391887        0.9055    108108        21        ATF7
#> 5             0.0391887        0.9055    108108        35       BRCA2
#> ...                 ...           ...       ...       ...         ...
#> 5214280       0.0420235       0.90550     44383      1366     ZSCAN21
#> 5214281       0.0420235       0.90550     44383      1376        ZXDB
#> 5214282       0.0057315       0.82698    105268       114       FOXA1
#> 5214283       0.0057315       0.82698    105268       128       GATA2
#> 5214284       0.0057315       0.82698    105268       823       SUMO2

4.6 Add Weights

While the pruneRegulon function provides statistics on the joint occurrence of TF-RE-TG, we would like to further estimate the strength of regulation. Biologically, this can be interpreted as the magnitude of gene expression changes induced by transcription factor activity. Epiregulon estimates the regulatory potential using one of the three measures: 1) correlation between TF and target gene expression, 2) mutual information between the TF and target gene expression and 3) Wilcoxon test statistics of target gene expression in cells jointly expressing all 3 elements vs cells that do not.

Two of the measures (correlation and Wilcoxon statistics) give both the magnitude and directionality of changes whereas mutual information is always positive. The correlation and mutual information statistics are computed on pseudobulks aggregated by user-supplied cluster labels, whereas the Wilcoxon method groups cells into two categories: 1) the active category of cells jointly expressing TF, RE and TG and 2) the inactive category consisting of the remaining cells.

We calculate cluster-specific weights if users supply cluster labels.


regulon.w <- addWeights(regulon = pruned.regulon,
                        expMatrix  = GeneExpressionMatrix,
                        exp_assay  = "normalizedCounts",
                        peakMatrix = PeakMatrix,
                        peak_assay = "counts",
                        clusters = GeneExpressionMatrix$Clusters,
                        method = "wilcox")

regulon.w

4.7 Annotate with TF motifs (Optional)

So far the gene regulatory network was constructed from TF ChIP-seq exclusively. Some users would prefer to further annotate regulatory elements with the presence of motifs. We provide an option to annotate peaks with motifs from the Cisbp database. If no motifs are present for this particular factor (as in the case of co-factors or chromatin modifiers), we return NA. If motifs are available for a factor and the RE contains a motif, we return 1. If motifs are available and the RE does not contain a motif, we return 0. The users can also provide their own motif annotation through the pwms argument.


regulon.w.motif <- addMotifScore(regulon = regulon.w,
                                 peaks = rowRanges(PeakMatrix),
                                 species = "human",
                                 genome = "hg38")
#> annotating peaks with motifs
#> see ?scMultiome and browseVignettes('scMultiome') for documentation
#> loading from cache
#> 
#> Attaching package: 'Biostrings'
#> The following object is masked from 'package:base':
#> 
#>     strsplit
#> 
#> Attaching package: 'rtracklayer'
#> The following object is masked from 'package:AnnotationHub':
#> 
#>     hubUrl

# if desired, set weight to 0 if no motif is found
regulon.w.motif$weight[regulon.w.motif$motif == 0] <- 0

regulon.w.motif
#> DataFrame with 12318 rows and 17 columns
#>         idxATAC         chr     start       end    idxRNA     target      corr
#>       <integer> <character> <integer> <integer> <integer>    <array>  <matrix>
#> 1            10        chr1    920987    921487        28      PERM1  0.602486
#> 2            14        chr1    932922    933422        28      PERM1 -0.408990
#> 3            14        chr1    932922    933422        20 AL645608.6 -0.486259
#> 4            23        chr1    966150    966650        33       AGRN -0.517957
#> 5            24        chr1    966766    967266        21 AL645608.2 -0.407948
#> ...         ...         ...       ...       ...       ...        ...       ...
#> 12314    122461        chr9 113174265 113174765     34825    SLC31A1  0.772916
#> 12315    123380        chr9 128635245 128635745     35034     SPTAN1  0.861189
#> 12316    124405        chrX  15737784  15738284     35402      AP1S2 -0.479022
#> 12317    124971        chrX  47217674  47218174     35588       UBA1 -0.475001
#> 12318    124978        chrX  47541057  47541557     35594      ZNF41 -0.456539
#>       p_val_peak_gene FDR_peak_gene  distance     idxTF          tf
#>              <matrix>      <matrix> <integer> <integer> <character>
#> 1          0.03013003      0.905500     59540       485          AR
#> 2          0.03068495      0.905500     47605       485          AR
#> 3          0.00973898      0.847422     28088       485          AR
#> 4          0.00592574      0.826980     53468       485          AR
#> 5          0.03118622      0.905500     55331       485          AR
#> ...               ...           ...       ...       ...         ...
#> 12314      0.00663767      0.826980     46777       723      NKX2-1
#> 12315      0.00217480      0.826980     82687       723      NKX2-1
#> 12316      0.01115328      0.864689    116469       723      NKX2-1
#> 12317      0.01181568      0.868919     26813       723      NKX2-1
#> 12318      0.01595116      0.898733     58472       723      NKX2-1
#>                            pval                    stats             qval
#>                        <matrix>                 <matrix>         <matrix>
#> 1     0.103465:0.88043479:1:... 2.651304:0.0226256:0:...        1:1:1:...
#> 2     0.559820:0.78785954:1:... 0.340016:0.0724101:0:...        1:1:1:...
#> 3     0.182226:0.50963499:1:... 1.779378:0.4348165:0:...        1:1:1:...
#> 4     0.238769:0.00883168:1:... 1.387848:6.8565532:0:...        1:1:1:...
#> 5     0.188330:0.66852363:1:... 1.730635:0.1833339:0:...        1:1:1:...
#> ...                         ...                      ...              ...
#> 12314       4.51202e-02:1:1:...          4.01414:0:0:... 1.000000:1:1:...
#> 12315       4.51202e-02:1:1:...          4.01414:0:0:... 1.000000:1:1:...
#> 12316       1.65845e-02:1:1:...          5.73982:0:0:... 1.000000:1:1:...
#> 12317       3.54310e-03:1:1:...          8.50429:0:0:... 1.000000:1:1:...
#> 12318       9.49637e-06:1:1:...         19.61014:0:0:... 0.674527:1:1:...
#>          weight     motif
#>        <matrix> <numeric>
#> 1     0:0:0:...         0
#> 2     0:0:0:...         0
#> 3     0:0:0:...         0
#> 4     0:0:0:...         0
#> 5     0:0:0:...         0
#> ...         ...       ...
#> 12314 0:0:0:...         0
#> 12315 0:0:0:...         0
#> 12316 0:0:0:...         0
#> 12317 0:0:0:...         0
#> 12318 0:0:0:...         0

4.8 Annotate with log fold changes (Optional)

It is sometimes helpful to filter the regulons based on gene expression changes between two conditions. The addLogFC function is a wrapper around scran::findMarkers and adds extra columns of log changes to the regulon DataFrame. The users can specify the reference condition in logFC_ref and conditions for comparison in logFC_condition. If these are not provided, log fold changes are calculated for every condition in the cluster labels against the rest of the conditions.

# create logcounts
GeneExpressionMatrix <- scuttle::logNormCounts(GeneExpressionMatrix)

# add log fold changes which are calculated by taking the difference of the log counts
regulon.w <- addLogFC(regulon = regulon.w,
                      clusters = GeneExpressionMatrix$hash_assignment,
                      expMatrix  = GeneExpressionMatrix,
                      assay.type  = "logcounts",
                      pval.type = "any",
                      logFC_condition = c("HTO10_GATA6_UTR", 
                                          "HTO2_GATA6_v2", 
                                          "HTO8_NKX2.1_UTR"),
                      logFC_ref = "HTO5_NeonG_v2")

regulon.w
#> DataFrame with 12318 rows and 25 columns
#>         idxATAC         chr     start       end    idxRNA     target      corr
#>       <integer> <character> <integer> <integer> <integer>    <array>  <matrix>
#> 1            10        chr1    920987    921487        28      PERM1  0.602486
#> 2            14        chr1    932922    933422        28      PERM1 -0.408990
#> 3            14        chr1    932922    933422        20 AL645608.6 -0.486259
#> 4            23        chr1    966150    966650        33       AGRN -0.517957
#> 5            24        chr1    966766    967266        21 AL645608.2 -0.407948
#> ...         ...         ...       ...       ...       ...        ...       ...
#> 12314    122461        chr9 113174265 113174765     34825    SLC31A1  0.772916
#> 12315    123380        chr9 128635245 128635745     35034     SPTAN1  0.861189
#> 12316    124405        chrX  15737784  15738284     35402      AP1S2 -0.479022
#> 12317    124971        chrX  47217674  47218174     35588       UBA1 -0.475001
#> 12318    124978        chrX  47541057  47541557     35594      ZNF41 -0.456539
#>       p_val_peak_gene FDR_peak_gene  distance     idxTF          tf
#>              <matrix>      <matrix> <integer> <integer> <character>
#> 1          0.03013003      0.905500     59540       485          AR
#> 2          0.03068495      0.905500     47605       485          AR
#> 3          0.00973898      0.847422     28088       485          AR
#> 4          0.00592574      0.826980     53468       485          AR
#> 5          0.03118622      0.905500     55331       485          AR
#> ...               ...           ...       ...       ...         ...
#> 12314      0.00663767      0.826980     46777       723      NKX2-1
#> 12315      0.00217480      0.826980     82687       723      NKX2-1
#> 12316      0.01115328      0.864689    116469       723      NKX2-1
#> 12317      0.01181568      0.868919     26813       723      NKX2-1
#> 12318      0.01595116      0.898733     58472       723      NKX2-1
#>                            pval                    stats             qval
#>                        <matrix>                 <matrix>         <matrix>
#> 1     0.103465:0.88043479:1:... 2.651304:0.0226256:0:...        1:1:1:...
#> 2     0.559820:0.78785954:1:... 0.340016:0.0724101:0:...        1:1:1:...
#> 3     0.182226:0.50963499:1:... 1.779378:0.4348165:0:...        1:1:1:...
#> 4     0.238769:0.00883168:1:... 1.387848:6.8565532:0:...        1:1:1:...
#> 5     0.188330:0.66852363:1:... 1.730635:0.1833339:0:...        1:1:1:...
#> ...                         ...                      ...              ...
#> 12314       4.51202e-02:1:1:...          4.01414:0:0:... 1.000000:1:1:...
#> 12315       4.51202e-02:1:1:...          4.01414:0:0:... 1.000000:1:1:...
#> 12316       1.65845e-02:1:1:...          5.73982:0:0:... 1.000000:1:1:...
#> 12317       3.54310e-03:1:1:...          8.50429:0:0:... 1.000000:1:1:...
#> 12318       9.49637e-06:1:1:...         19.61014:0:0:... 0.674527:1:1:...
#>                             weight HTO10_GATA6_UTR.vs.HTO5_NeonG_v2.FDR
#>                           <matrix>                            <numeric>
#> 1     0.02638820:-0.00723681:0:...                            1.0000000
#> 2     0.00967161:-0.01328838:0:...                            1.0000000
#> 3     0.03590731:-0.03273486:0:...                            0.3913742
#> 4     0.02697100: 0.15302114:0:...                            0.0333213
#> 5     0.02780082:-0.02106160:0:...                            0.7803240
#> ...                            ...                                  ...
#> 12314            0.0333575:0:0:...                            0.0771375
#> 12315            0.0277381:0:0:...                            0.8910514
#> 12316            0.0397841:0:0:...                            1.0000000
#> 12317            0.0449399:0:0:...                            1.0000000
#> 12318            0.0696895:0:0:...                            1.0000000
#>       HTO10_GATA6_UTR.vs.HTO5_NeonG_v2.p.value
#>                                      <numeric>
#> 1                                    -0.306482
#> 2                                    -0.306482
#> 3                                    -3.630700
#> 4                                    -6.908244
#> 5                                    -2.358197
#> ...                                        ...
#> 12314                               -5.8525011
#> 12315                               -2.0888086
#> 12316                               -0.8173609
#> 12317                               -0.0600713
#> 12318                               -1.0486948
#>       HTO10_GATA6_UTR.vs.HTO5_NeonG_v2.logFC HTO2_GATA6_v2.vs.HTO5_NeonG_v2.FDR
#>                                    <numeric>                          <numeric>
#> 1                                 0.00267018                         0.95051180
#> 2                                 0.00267018                         0.95051180
#> 3                                -0.02908575                         1.00000000
#> 4                                 0.06228930                         0.00813329
#> 5                                -0.01721007                         1.00000000
#> ...                                      ...                                ...
#> 12314                            -0.09733426                          0.6730897
#> 12315                             0.05722383                          0.0780338
#> 12316                            -0.01490592                          0.9505118
#> 12317                            -0.00210678                          1.0000000
#> 12318                             0.01652101                          1.0000000
#>       HTO2_GATA6_v2.vs.HTO5_NeonG_v2.p.value
#>                                    <numeric>
#> 1                                  -1.145473
#> 2                                  -1.145473
#> 3                                  -0.613957
#> 4                                  -8.744274
#> 5                                  -0.916768
#> ...                                      ...
#> 12314                              -2.795089
#> 12315                              -5.984106
#> 12316                              -1.533481
#> 12317                              -0.668397
#> 12318                              -0.434791
#>       HTO2_GATA6_v2.vs.HTO5_NeonG_v2.logFC HTO8_NKX2.1_UTR.vs.HTO5_NeonG_v2.FDR
#>                                  <numeric>                            <numeric>
#> 1                              -0.00410095                             1.000000
#> 2                              -0.00410095                             1.000000
#> 3                              -0.00903561                             1.000000
#> 4                               0.06470266                             0.968083
#> 5                              -0.00918756                             1.000000
#> ...                                    ...                                  ...
#> 12314                          -0.05908672                             0.991399
#> 12315                           0.11104599                             0.992412
#> 12316                          -0.02184636                             0.978457
#> 12317                          -0.01666954                             1.000000
#> 12318                           0.00696425                             0.978457
#>       HTO8_NKX2.1_UTR.vs.HTO5_NeonG_v2.p.value
#>                                      <numeric>
#> 1                                    -0.125877
#> 2                                    -0.125877
#> 3                                    -0.318033
#> 4                                    -2.479073
#> 5                                    -0.582022
#> ...                                        ...
#> 12314                                -1.124680
#> 12315                                -1.122331
#> 12316                                -1.855533
#> 12317                                -0.689542
#> 12318                                -1.356691
#>       HTO8_NKX2.1_UTR.vs.HTO5_NeonG_v2.logFC
#>                                    <numeric>
#> 1                               -0.000708047
#> 2                               -0.000708047
#> 3                                0.005164618
#> 4                                0.024761627
#> 5                               -0.006427500
#> ...                                      ...
#> 12314                              0.0332807
#> 12315                              0.0314664
#> 12316                             -0.0243966
#> 12317                             -0.0168862
#> 12318                              0.0170046

4.9 Calculate TF activity

Finally, the activities for a specific TF in each cell are computed by averaging expressions of target genes linked to the TF weighted by the test statistics of choice, chosen from either correlation, mutual information or the Wilcoxon test statistics. \[y=\frac{1}{n}\sum_{i=1}^{n} x_i * weights_i\] where \(y\) is the activity of a TF for a cell, \(n\) is the total number of targets for a TF, \(x_i\) is the log count expression of target \(i\) where \(i\) in {1,2,…,n} and \(weights_i\) is the weight of TF - target \(i\)

score.combine <- calculateActivity(expMatrix = GeneExpressionMatrix,
                                   regulon = regulon.w, 
                                   mode = "weight", 
                                   method = "weightedMean", 
                                   exp_assay = "normalizedCounts",
                                   normalize = FALSE)
#> Warning in calculateActivity(expMatrix = GeneExpressionMatrix, regulon =
#> regulon.w, : Argument 'method' to calculateActivity was deprecated as of
#> epiregulon version 2.0.0
#> calculating TF activity from regulon using weightedMean
#> Warning in calculateActivity(expMatrix = GeneExpressionMatrix, regulon = regulon.w, : The weight column contains multiple subcolumns but no cluster information was provided.
#>           Using first column to compute activity...
#> aggregating regulons...
#> creating weight matrix...
#> calculating activity scores...
#> normalize by the number of targets...

score.combine[1:5,1:5]
#> 5 x 5 sparse Matrix of class "dgCMatrix"
#>        reprogram#TTAGGAACAAGGTACG-1 reprogram#GAGCGGTCAACCTGGT-1
#> AR                       0.01775226                   0.02035383
#> FOXA1                    0.01671665                   0.01742678
#> FOXA2                    0.01589522                   0.02287218
#> GATA6                    0.01531109                   0.05477733
#> NKX2-1                   0.03161050                   0.01936065
#>        reprogram#TTATAGCCACCCTCAC-1 reprogram#TGGTGATTCCTGTTCA-1
#> AR                      0.010839366                  0.010916741
#> FOXA1                   0.010880353                  0.011025560
#> FOXA2                   0.012400559                  0.013676874
#> GATA6                   0.008777784                  0.007779859
#> NKX2-1                  0.013159392                  0.011983253
#>        reprogram#TCGGTTCTCACTAGGT-1
#> AR                       0.01810095
#> FOXA1                    0.01644120
#> FOXA2                    0.01640028
#> GATA6                    0.01299646
#> NKX2-1                   0.02691256

5 Session Info

sessionInfo()
#> R Under development (unstable) (2025-10-20 r88955)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.3 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.23-bioc/R/lib/libRblas.so 
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0  LAPACK version 3.12.0
#> 
#> 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       
#> 
#> time zone: America/New_York
#> tzcode source: system (glibc)
#> 
#> attached base packages:
#> [1] stats4    stats     graphics  grDevices utils     datasets  methods  
#> [8] base     
#> 
#> other attached packages:
#>  [1] BSgenome.Hsapiens.UCSC.hg38_1.4.5 BSgenome_1.79.1                  
#>  [3] rtracklayer_1.71.0                BiocIO_1.21.0                    
#>  [5] Biostrings_2.79.2                 XVector_0.51.0                   
#>  [7] GenomeInfoDb_1.47.0               scMultiome_1.11.0                
#>  [9] MultiAssayExperiment_1.37.1       ExperimentHub_3.1.0              
#> [11] AnnotationHub_4.1.0               BiocFileCache_3.1.0              
#> [13] dbplyr_2.5.1                      epiregulon_2.1.2                 
#> [15] SingleCellExperiment_1.33.0       SummarizedExperiment_1.41.0      
#> [17] Biobase_2.71.0                    GenomicRanges_1.63.0             
#> [19] Seqinfo_1.1.0                     IRanges_2.45.0                   
#> [21] S4Vectors_0.49.0                  BiocGenerics_0.57.0              
#> [23] generics_0.1.4                    MatrixGenerics_1.23.0            
#> [25] matrixStats_1.5.0                 BiocStyle_2.39.0                 
#> 
#> loaded via a namespace (and not attached):
#>   [1] RColorBrewer_1.1-3          jsonlite_2.0.0             
#>   [3] magrittr_2.0.4              ggbeeswarm_0.7.2           
#>   [5] magick_2.9.0                farver_2.1.2               
#>   [7] rmarkdown_2.30              vctrs_0.6.5                
#>   [9] memoise_2.0.1               Rsamtools_2.27.0           
#>  [11] RCurl_1.98-1.17             tinytex_0.57               
#>  [13] htmltools_0.5.8.1           S4Arrays_1.11.0            
#>  [15] curl_7.0.0                  BiocNeighbors_2.5.0        
#>  [17] Rhdf5lib_1.33.0             SparseArray_1.11.1         
#>  [19] rhdf5_2.55.8                sass_0.4.10                
#>  [21] bslib_0.9.0                 httr2_1.2.1                
#>  [23] cachem_1.1.0                GenomicAlignments_1.47.0   
#>  [25] igraph_2.2.1                lifecycle_1.0.4            
#>  [27] pkgconfig_2.0.3             rsvd_1.0.5                 
#>  [29] Matrix_1.7-4                R6_2.6.1                   
#>  [31] fastmap_1.2.0               digest_0.6.38              
#>  [33] TFMPvalue_0.0.9             AnnotationDbi_1.73.0       
#>  [35] scater_1.39.0               dqrng_0.4.1                
#>  [37] irlba_2.3.5.1               RSQLite_2.4.4              
#>  [39] beachmat_2.27.0             seqLogo_1.77.0             
#>  [41] filelock_1.0.3              labeling_0.4.3             
#>  [43] httr_1.4.7                  abind_1.4-8                
#>  [45] compiler_4.6.0              bit64_4.6.0-1              
#>  [47] withr_3.0.2                 S7_0.2.0                   
#>  [49] backports_1.5.0             BiocParallel_1.45.0        
#>  [51] viridis_0.6.5               DBI_1.2.3                  
#>  [53] HDF5Array_1.39.0            rappdirs_0.3.3             
#>  [55] DelayedArray_0.37.0         bluster_1.21.0             
#>  [57] rjson_0.2.23                gtools_3.9.5               
#>  [59] caTools_1.18.3              tools_4.6.0                
#>  [61] vipor_0.4.7                 beeswarm_0.4.0             
#>  [63] glue_1.8.0                  h5mread_1.3.0              
#>  [65] restfulr_0.0.16             rhdf5filters_1.23.0        
#>  [67] grid_4.6.0                  checkmate_2.3.3            
#>  [69] cluster_2.1.8.1             TFBSTools_1.49.0           
#>  [71] gtable_0.3.6                metapod_1.19.0             
#>  [73] BiocSingular_1.27.0         ScaledMatrix_1.19.0        
#>  [75] ggrepel_0.9.6               BiocVersion_3.23.1         
#>  [77] pillar_1.11.1               motifmatchr_1.33.0         
#>  [79] limma_3.67.0                dplyr_1.1.4                
#>  [81] lattice_0.22-7              bit_4.6.0                  
#>  [83] DirichletMultinomial_1.53.0 tidyselect_1.2.1           
#>  [85] locfit_1.5-9.12             scuttle_1.21.0             
#>  [87] knitr_1.50                  gridExtra_2.3              
#>  [89] scrapper_1.5.2              bookdown_0.45              
#>  [91] edgeR_4.9.0                 xfun_0.54                  
#>  [93] statmod_1.5.1               UCSC.utils_1.7.0           
#>  [95] yaml_2.3.10                 evaluate_1.0.5             
#>  [97] codetools_0.2-20            cigarillo_1.1.0            
#>  [99] tibble_3.3.0                BiocManager_1.30.26        
#> [101] cli_3.6.5                   jquerylib_0.1.4            
#> [103] dichromat_2.0-0.1           Rcpp_1.1.0                 
#> [105] png_0.1-8                   XML_3.99-0.20              
#> [107] parallel_4.6.0              ggplot2_4.0.0              
#> [109] blob_1.2.4                  scran_1.39.0               
#> [111] bitops_1.0-9                pwalign_1.7.0              
#> [113] viridisLite_0.4.2           scales_1.4.0               
#> [115] purrr_1.2.0                 crayon_1.5.3               
#> [117] rlang_1.1.6                 cowplot_1.2.0              
#> [119] KEGGREST_1.51.0