Contents

1 Introduction

epiregulon.extra is a companion package to epiregulon and provides functions to visualize transcription factor activity and network. It also provides statistical tests to identify transcription factors with differential activity and network topology. This tutorial continues from the reprogram-seq example included in epiregulon. We will use the gene regulatory networks constructed by epiregulon and continue with visualization and network analysis.

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.extra")

3 Data preparation

To continue with the visualization functions, we will first need the gene expression matrix. We can download the data from scMultiome.

# load the MAE object
library(scMultiome)
## Loading required package: AnnotationHub
## Loading required package: BiocGenerics
## 
## Attaching package: 'BiocGenerics'
## The following objects are masked from 'package:stats':
## 
##     IQR, mad, sd, var, xtabs
## The following objects are masked from 'package:base':
## 
##     Filter, Find, Map, Position, Reduce, anyDuplicated, aperm, append,
##     as.data.frame, basename, cbind, colnames, dirname, do.call,
##     duplicated, eval, evalq, get, grep, grepl, intersect, is.unsorted,
##     lapply, mapply, match, mget, order, paste, pmax, pmax.int, pmin,
##     pmin.int, rank, rbind, rownames, sapply, setdiff, table, tapply,
##     union, unique, unsplit, which.max, which.min
## Loading required package: BiocFileCache
## Loading required package: dbplyr
## Loading required package: ExperimentHub
## Loading required package: MultiAssayExperiment
## Loading required package: SummarizedExperiment
## Loading required package: MatrixGenerics
## Loading required package: matrixStats
## 
## Attaching package: 'MatrixGenerics'
## The following objects are masked from 'package:matrixStats':
## 
##     colAlls, colAnyNAs, colAnys, colAvgsPerRowSet, colCollapse,
##     colCounts, colCummaxs, colCummins, colCumprods, colCumsums,
##     colDiffs, colIQRDiffs, colIQRs, colLogSumExps, colMadDiffs,
##     colMads, colMaxs, colMeans2, colMedians, colMins, colOrderStats,
##     colProds, colQuantiles, colRanges, colRanks, colSdDiffs, colSds,
##     colSums2, colTabulates, colVarDiffs, colVars, colWeightedMads,
##     colWeightedMeans, colWeightedMedians, colWeightedSds,
##     colWeightedVars, rowAlls, rowAnyNAs, rowAnys, rowAvgsPerColSet,
##     rowCollapse, rowCounts, rowCummaxs, rowCummins, rowCumprods,
##     rowCumsums, rowDiffs, rowIQRDiffs, rowIQRs, rowLogSumExps,
##     rowMadDiffs, rowMads, rowMaxs, rowMeans2, rowMedians, rowMins,
##     rowOrderStats, rowProds, rowQuantiles, rowRanges, rowRanks,
##     rowSdDiffs, rowSds, rowSums2, rowTabulates, rowVarDiffs, rowVars,
##     rowWeightedMads, rowWeightedMeans, rowWeightedMedians,
##     rowWeightedSds, rowWeightedVars
## Loading required package: GenomicRanges
## Loading required package: stats4
## Loading required package: S4Vectors
## 
## Attaching package: 'S4Vectors'
## The following object is masked from 'package:utils':
## 
##     findMatches
## The following objects are masked from 'package:base':
## 
##     I, expand.grid, unname
## Loading required package: IRanges
## Loading required package: GenomeInfoDb
## Loading required package: Biobase
## Welcome to Bioconductor
## 
##     Vignettes contain introductory material; view with
##     'browseVignettes()'. To cite Bioconductor, see
##     'citation("Biobase")', and for packages 'citation("pkgname")'.
## 
## Attaching package: 'Biobase'
## The following object is masked from 'package:MatrixGenerics':
## 
##     rowMedians
## The following objects are masked from 'package:matrixStats':
## 
##     anyMissing, rowMedians
## The following object is masked from 'package:ExperimentHub':
## 
##     cache
## The following object is masked from 'package:AnnotationHub':
## 
##     cache
## Loading required package: SingleCellExperiment
mae <- scMultiome::reprogramSeq()
## see ?scMultiome and browseVignettes('scMultiome') for documentation
## loading from cache
# expression matrix
GeneExpressionMatrix <- mae[["GeneExpressionMatrix"]]
rownames(GeneExpressionMatrix) <- rowData(GeneExpressionMatrix)$name

reducedDim(GeneExpressionMatrix, "UMAP_Combined") <- 
  reducedDim(mae[["TileMatrix500"]], "UMAP_Combined")

# arrange hash_assignment
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"
  )
)

4 Calculate TF activity

Next we load the regulon object previously calculated by epiregulon. Here we just use the pruned version of regulon object in which only relevant columns are kept. Using epiregulon, we can calculate activity of transcription factors included in the regulon object.

library(epiregulon)
library(epiregulon.extra)

data(regulon)

score.combine <- calculateActivity(
  expMatrix = GeneExpressionMatrix,
  regulon = regulon,
  mode = "weight",
  method = "weightedMean",
  exp_assay = "normalizedCounts",
  normalize = FALSE
)
## calculating TF activity from regulon using weightedmean
## Warning in calculateActivity(expMatrix = GeneExpressionMatrix, regulon =
## regulon, : 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...

5 Perform differential activity

We perform a differential analysis of transcription factor activity across groups of cells. This function is a wrapper around findMarkers from scran.

library(epiregulon.extra)
markers <- findDifferentialActivity(
  activity_matrix = score.combine,
  clusters = GeneExpressionMatrix$hash_assignment,
  pval.type = "some",
  direction = "up",
  test.type = "t"
)

We can specify the top transcription factors of interest

markers.sig <- getSigGenes(markers, topgenes = 5)
## Using a logFC cutoff of 0 for class HTO10_GATA6_UTR for direction equal to any
## Using a logFC cutoff of 0 for class HTO2_GATA6_v2 for direction equal to any
## Using a logFC cutoff of 0 for class HTO8_NKX2.1_UTR for direction equal to any
## Using a logFC cutoff of 0 for class HTO3_NKX2.1_v2 for direction equal to any
## Using a logFC cutoff of 0 for class HTO1_FOXA2_v2 for direction equal to any
## Using a logFC cutoff of 0 for class HTO4_mFOXA1_v2 for direction equal to any
## Using a logFC cutoff of 0 for class HTO6_hFOXA1_UTR for direction equal to any
## Using a logFC cutoff of 0 for class HTO5_NeonG_v2 for direction equal to any

6 Visualize the results

We visualize the known differential transcription factors by bubble plot

plotBubble(
  activity_matrix = score.combine,
  tf = c("NKX2-1", "GATA6", "FOXA1", "FOXA2"),
  clusters = GeneExpressionMatrix$hash_assignment
)

Then visualize the most differential transcription factors by clusters

plotBubble(
  activity_matrix = score.combine,
  tf = markers.sig$tf,
  clusters = GeneExpressionMatrix$hash_assignment
)

Visualize the known differential transcription factors by violin plot.

plotActivityViolin(
  activity_matrix = score.combine,
  tf = c("NKX2-1", "GATA6", "FOXA1", "FOXA2", "AR"),
  clusters = GeneExpressionMatrix$hash_assignment
)

Visualize the known differential transcription factors by UMAP

options(ggrastr.default.dpi=100)

ActivityPlot <- plotActivityDim(
  sce = GeneExpressionMatrix,
  activity_matrix = score.combine,
  tf = c("NKX2-1", "GATA6", "FOXA1", "FOXA2", "AR"),
  dimtype = "UMAP_Combined",
  label = "Clusters",
  point_size = 0.3,
  ncol = 3,
  rasterise = TRUE
)


ActivityPlot

In contrast, the gene expression of the transcription factors is very sparse

options(ggrastr.default.dpi=100)

ActivityPlot <- plotActivityDim(
  sce = GeneExpressionMatrix,
  activity_matrix = counts(GeneExpressionMatrix),
  tf = c("NKX2-1", "GATA6", "FOXA1", "FOXA2", "AR"),
  dimtype = "UMAP_Combined",
  label = "Clusters",
  point_size = 0.3,
  ncol = 3,
  limit = c(0, 2),
  colors = c("grey", "blue"),
  legend.label = "GEX",
  rasterise = TRUE
)
ActivityPlot

Visualize the activity of the selected transcription factors by heat map. Due to the maximum size limit for this vignette, the output is not shown here.

plotHeatmapRegulon(
  sce = GeneExpressionMatrix,
  tfs = c("GATA6", "NKX2-1"),
  regulon = regulon,
  regulon_cutoff = 0,
  downsample = 1000,
  cell_attributes = "hash_assignment",
  col_gap = "hash_assignment",
  exprs_values = "normalizedCounts",
  name = "regulon heatmap",
  column_title_rot = 45,
  raster_quality=4
)


plotHeatmapActivity(
  activity = score.combine,
  sce = GeneExpressionMatrix,
  tfs = c("GATA6", "NKX2-1", "FOXA1", "FOXA2"),
  downsample = 5000,
  cell_attributes = "hash_assignment",
  col_gap = "hash_assignment",
  name = "Activity",
  column_title_rot = 45,
  raster_quality=3
)

7 Geneset enrichment

Sometimes we are interested to know what pathways are enriched in the regulon of a particular TF. We can perform geneset enrichment using the enricher function from clusterProfiler.

# retrieve genesets
H <- EnrichmentBrowser::getGenesets(
  org = "hsa",
  db = "msigdb",
  cat = "H",
  gene.id.type = "SYMBOL"
)
C2 <- EnrichmentBrowser::getGenesets(
  org = "hsa",
  db = "msigdb",
  cat = "C2",
  gene.id.type = "SYMBOL"
)
C6 <- EnrichmentBrowser::getGenesets(
  org = "hsa",
  db = "msigdb",
  cat = "C6",
  gene.id.type = "SYMBOL"
)

# combine genesets and convert genesets to be compatible with enricher
gs <- c(H, C2, C6)
gs.list <- do.call(rbind, lapply(names(gs), function(x) {
  data.frame(gs = x, genes = gs[[x]])
}))

enrichresults <- regulonEnrich(
  TF = c("GATA6", "NKX2-1"),
  regulon = regulon,
  weight = "weight",
  weight_cutoff = 0.1,
  genesets = gs.list
)
## GATA6
## 
## NKX2-1
# plot results
enrichPlot(results = enrichresults)

8 Network analysis

We can visualize the genesets as a network

plotGseaNetwork(
  tf = names(enrichresults),
  enrichresults = enrichresults,
  p.adj_cutoff = 0.1,
  ntop_pathways = 10
)

9 Differential networks

We are interested in understanding the differential networks between two conditions and determining which transcription factors account for the differences in the topology of networks. The pruned regulons with cluster-specific test statistics computed by pruneRegulon can be used to generate cluster-specific networks based on user-defined cutoffs and to visualize differential networks for transcription factors of interest. In this dataset, the GATA6 gene was only expressed in cluster 1 (C1) and NKX2-1 was only expressed in cluster 3 (C3). If we visualize the target genes of GATA6, we can see that C1 has many more target genes of GATA6 compared to C5, a cluster that does not express GATA6. Similarly, NKX2-1 target genes are confined to C3 which is the only cluster that exogenously expresses NKX2-1.

plotDiffNetwork(regulon,
  cutoff = 0.1,
  tf = c("GATA6"),
  weight = "weight",
  clusters = c("C1", "C5"),
  layout = "stress"
)
## Building graph using weight as edge weights

plotDiffNetwork(regulon,
  cutoff = 0.1,
  tf = c("NKX2-1"),
  weight = "weight",
  clusters = c("C3", "C5"),
  layout = "stress"
)
## Building graph using weight as edge weights

We can also visualize how transcription factors relate to other transcription factors in C1 cluster.

selected <- which(regulon$weight[, "C1"] > 0 &
  regulon$tf %in% c("GATA6", "FOXA1", "AR"))

C1_network <- buildGraph(regulon[selected, ],
  weights = "weight",
  cluster = "C1"
  )
## Building graph using weight as edge weights
plotEpiregulonNetwork(C1_network,
  layout = "sugiyama",
  tfs_to_highlight = c("GATA6", "FOXA1", "AR")) + 
  ggplot2::ggtitle("C1")

To systematically examine the differential network topology between two clusters, we perform an edge subtraction between two graphs, using weights computed by pruneRegulon. We then calculate the degree centrality of the weighted differential graphs and if desired, normalize the differential centrality against the total number of edges. The default normalization function is sqrt as it preserves both the difference in the number of edges (but scaled by sqrt) and the differences in the weights. If the user only wants to examine the differences in the averaged weights, the FUN argument can be changed to identity. Finally, we rank the transcription factors by (normalized) differential centrality. For demonstration purpose, the regulon list is truncated but the full list would contain all the transcription factors.

# rank by differential centrality
C1_network <- buildGraph(regulon, weights = "weight", cluster = "C1")
## Building graph using weight as edge weights
C5_network <- buildGraph(regulon, weights = "weight", cluster = "C5")
## Building graph using weight as edge weights
diff_graph <- buildDiffGraph(C1_network, C5_network, abs_diff = FALSE)
diff_graph <- addCentrality(diff_graph)
diff_graph <- normalizeCentrality(diff_graph)
rank_table <- rankTfs(diff_graph)

library(ggplot2)
ggplot(rank_table, aes(x = rank, y = centrality)) +
  geom_point() +
  ggrepel::geom_text_repel(data = rank_table, aes(label = tf)) +
  theme_classic()

10 Session Info

sessionInfo()
## R version 4.4.1 (2024-06-14)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 22.04.4 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.19-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.10.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] ggplot2_3.5.1               org.Hs.eg.db_3.19.1        
##  [3] AnnotationDbi_1.66.0        msigdbr_7.5.1              
##  [5] epiregulon.extra_1.0.1      epiregulon_1.0.1           
##  [7] scMultiome_1.4.0            SingleCellExperiment_1.26.0
##  [9] MultiAssayExperiment_1.30.3 SummarizedExperiment_1.34.0
## [11] Biobase_2.64.0              GenomicRanges_1.56.1       
## [13] GenomeInfoDb_1.40.1         IRanges_2.38.1             
## [15] S4Vectors_0.42.1            MatrixGenerics_1.16.0      
## [17] matrixStats_1.3.0           ExperimentHub_2.12.0       
## [19] AnnotationHub_3.12.0        BiocFileCache_2.12.0       
## [21] dbplyr_2.5.0                BiocGenerics_0.50.0        
## [23] BiocStyle_2.32.1           
## 
## loaded via a namespace (and not attached):
##   [1] splines_4.4.1             ggplotify_0.1.2          
##   [3] bitops_1.0-7              filelock_1.0.3           
##   [5] tibble_3.2.1              polyclip_1.10-6          
##   [7] graph_1.82.0              XML_3.99-0.17            
##   [9] lifecycle_1.0.4           edgeR_4.2.1              
##  [11] lattice_0.22-6            MASS_7.3-61              
##  [13] backports_1.5.0           magrittr_2.0.3           
##  [15] limma_3.60.3              sass_0.4.9               
##  [17] rmarkdown_2.27            jquerylib_0.1.4          
##  [19] yaml_2.3.9                metapod_1.12.0           
##  [21] RColorBrewer_1.1-3        cowplot_1.1.3            
##  [23] DBI_1.2.3                 abind_1.4-5              
##  [25] zlibbioc_1.50.0           purrr_1.0.2              
##  [27] ggraph_2.2.1              RCurl_1.98-1.16          
##  [29] yulab.utils_0.1.4         tweenr_2.0.3             
##  [31] rappdirs_0.3.3            GenomeInfoDbData_1.2.12  
##  [33] enrichplot_1.24.0         ggrepel_0.9.5            
##  [35] irlba_2.3.5.1             tidytree_0.4.6           
##  [37] annotate_1.82.0           dqrng_0.4.1              
##  [39] DelayedMatrixStats_1.26.0 codetools_0.2-20         
##  [41] DelayedArray_0.30.1       DOSE_3.30.1              
##  [43] scuttle_1.14.0            ggforce_0.4.2            
##  [45] tidyselect_1.2.1          aplot_0.2.3              
##  [47] UCSC.utils_1.0.0          farver_2.1.2             
##  [49] ScaledMatrix_1.12.0       viridis_0.6.5            
##  [51] jsonlite_1.8.8            BiocNeighbors_1.22.0     
##  [53] tidygraph_1.3.1           scater_1.32.0            
##  [55] tools_4.4.1               treeio_1.28.0            
##  [57] Rcpp_1.0.12               glue_1.7.0               
##  [59] gridExtra_2.3             SparseArray_1.4.8        
##  [61] xfun_0.45                 qvalue_2.36.0            
##  [63] dplyr_1.1.4               HDF5Array_1.32.0         
##  [65] withr_3.0.0               BiocManager_1.30.23      
##  [67] fastmap_1.2.0             rhdf5filters_1.16.0      
##  [69] bluster_1.14.0            fansi_1.0.6              
##  [71] digest_0.6.36             rsvd_1.0.5               
##  [73] gridGraphics_0.5-1        R6_2.5.1                 
##  [75] mime_0.12                 colorspace_2.1-0         
##  [77] GO.db_3.19.1              Cairo_1.6-2              
##  [79] RSQLite_2.3.7             utf8_1.2.4               
##  [81] tidyr_1.3.1               generics_0.1.3           
##  [83] data.table_1.15.4         graphlayouts_1.1.1       
##  [85] httr_1.4.7                S4Arrays_1.4.1           
##  [87] scatterpie_0.2.3          pkgconfig_2.0.3          
##  [89] gtable_0.3.5              blob_1.2.4               
##  [91] XVector_0.44.0            shadowtext_0.1.3         
##  [93] clusterProfiler_4.12.0    htmltools_0.5.8.1        
##  [95] fgsea_1.30.0              bookdown_0.40            
##  [97] GSEABase_1.66.0           scales_1.3.0             
##  [99] png_0.1-8                 ggfun_0.1.5              
## [101] scran_1.32.0              knitr_1.48               
## [103] reshape2_1.4.4            nlme_3.1-165             
## [105] checkmate_2.3.1           curl_5.2.1               
## [107] cachem_1.1.0              rhdf5_2.48.0             
## [109] stringr_1.5.1             BiocVersion_3.19.1       
## [111] HDO.db_0.99.1             parallel_4.4.1           
## [113] vipor_0.4.7               ggrastr_1.0.2            
## [115] pillar_1.9.0              grid_4.4.1               
## [117] vctrs_0.6.5               BiocSingular_1.20.0      
## [119] beachmat_2.20.0           xtable_1.8-4             
## [121] cluster_2.1.6             beeswarm_0.4.0           
## [123] Rgraphviz_2.48.0          evaluate_0.24.0          
## [125] KEGGgraph_1.64.0          tinytex_0.51             
## [127] magick_2.8.4              cli_3.6.3                
## [129] locfit_1.5-9.10           compiler_4.4.1           
## [131] rlang_1.1.4               crayon_1.5.3             
## [133] labeling_0.4.3            fs_1.6.4                 
## [135] plyr_1.8.9                ggbeeswarm_0.7.2         
## [137] stringi_1.8.4             viridisLite_0.4.2        
## [139] BiocParallel_1.38.0       babelgene_22.9           
## [141] munsell_0.5.1             Biostrings_2.72.1        
## [143] lazyeval_0.2.2            GOSemSim_2.30.0          
## [145] Matrix_1.7-0              patchwork_1.2.0          
## [147] sparseMatrixStats_1.16.0  bit64_4.0.5              
## [149] Rhdf5lib_1.26.0           KEGGREST_1.44.1          
## [151] statmod_1.5.0             highr_0.11               
## [153] igraph_2.0.3              memoise_2.0.1            
## [155] bslib_0.7.0               ggtree_3.12.0            
## [157] fastmatch_1.1-4           bit_4.0.5                
## [159] EnrichmentBrowser_2.34.1  gson_0.1.0               
## [161] ape_5.8