Bulk RNA-seq yield many molecular readouts that are hard to interpret by themselves. One way of summarizing this information is by inferring pathway activities from prior knowledge.

In this notebook we showcase how to use decoupleR for pathway activity inference with a bulk RNA-seq data-set where the transcription factor FOXA2 was knocked out in pancreatic cancer cell lines.

The data consists of 3 Wild Type (WT) samples and 3 Knock Outs (KO). They are freely available in GEO.

1 Loading packages

First, we need to load the relevant packages:

## We load the required packages
library(decoupleR)
library(dplyr)
library(tibble)
library(tidyr)
library(ggplot2)
library(pheatmap)
library(ggrepel)

2 Loading the data-set

Here we used an already processed bulk RNA-seq data-set. We provide the normalized log-transformed counts, the experimental design meta-data and the Differential Expressed Genes (DEGs) obtained using limma. We can open the data like this:

inputs_dir <- system.file("extdata", package = "decoupleR")
data <- readRDS(file.path(inputs_dir, "bk_data.rds"))

From data we can extract the mentioned information. Here we see the normalized log-transformed counts:

# Remove NAs and set row names
counts <- data$counts %>%
  dplyr::mutate_if(~ any(is.na(.x)), ~ if_else(is.na(.x),0,.x)) %>% 
  column_to_rownames(var = "gene") %>% 
  as.matrix()
head(counts)
#>          PANC1.WT.Rep1 PANC1.WT.Rep2 PANC1.WT.Rep3 PANC1.FOXA2KO.Rep1 PANC1.FOXA2KO.Rep2 PANC1.FOXA2KO.Rep3
#> NOC2L        10.052588     11.949123     12.057774          12.312291          12.139918          11.494205
#> PLEKHN1       7.535115      8.125993      8.714880           8.048196           8.290154           8.621239
#> PERM1         6.281242      6.424582      6.589668           6.293285           6.486136           6.775344
#> ISG15        10.938252     11.469081     11.425415          11.549986          11.371464          11.178157
#> AGRN          6.956335      7.196108      7.522550           7.061549           7.485534           7.071555
#> C1orf159      9.546224      9.788721      9.794589           9.850830           9.988069           9.965357

The design meta-data:

design <- data$design
design
#> # A tibble: 6 × 2
#>   sample             condition    
#>   <chr>              <chr>        
#> 1 PANC1.WT.Rep1      PANC1.WT     
#> 2 PANC1.WT.Rep2      PANC1.WT     
#> 3 PANC1.WT.Rep3      PANC1.WT     
#> 4 PANC1.FOXA2KO.Rep1 PANC1.FOXA2KO
#> 5 PANC1.FOXA2KO.Rep2 PANC1.FOXA2KO
#> 6 PANC1.FOXA2KO.Rep3 PANC1.FOXA2KO

And the results of limma, of which we are interested in extracting the obtained t-value from the contrast:

# Extract t-values per gene
deg <- data$limma_ttop %>%
    select(ID, t) %>% 
    filter(!is.na(t)) %>% 
    column_to_rownames(var = "ID") %>%
    as.matrix()
head(deg)
#>                  t
#> RHBDL2  -12.810588
#> PLEKHH2 -10.794453
#> HEG1     -9.788112
#> CLU      -9.761618
#> FHL1      8.950191
#> RBP4     -8.529074

3 PROGENy model

PROGENy is a comprehensive resource containing a curated collection of pathways and their target genes, with weights for each interaction. For this example we will use the human weights (mouse is also available) and we will use the top 100 responsive genes ranked by p-value. We can use decoupleR to retrieve it from OmniPath:

net <- get_progeny(organism = 'human', top = 100)
net
#> # A tibble: 1,400 × 4
#>    source   target  weight  p_value
#>    <chr>    <chr>    <dbl>    <dbl>
#>  1 Androgen TMPRSS2  11.5  2.38e-47
#>  2 Androgen NKX3-1   10.6  2.21e-44
#>  3 Androgen MBOAT2   10.5  4.63e-44
#>  4 Androgen KLK2     10.2  1.94e-40
#>  5 Androgen SARG     11.4  2.79e-40
#>  6 Androgen SLC38A4   7.36 1.25e-39
#>  7 Androgen MTMR9     6.13 2.53e-38
#>  8 Androgen ZBTB16   10.6  1.57e-36
#>  9 Androgen KCNN2     9.47 7.71e-36
#> 10 Androgen OPRK1    -5.63 1.11e-35
#> # … with 1,390 more rows

4 Activity inference with Weighted Mean

To infer activities we will run the Weighted Mean method (wmean). It infers regulator activities by first multiplying each target feature by its associated weight which then are summed to an enrichment score wmean. Furthermore, permutations of random target features can be performed to obtain a null distribution that can be used to compute a z-score norm_wmean, or a corrected estimate corr_wmean by multiplying wmean by the minus log10 of the obtained empirical p-value.

In this example we use wmean but we could have used any other. To see what methods are available use show_methods().

To run decoupleR methods, we need an input matrix (mat), an input prior knowledge network/resource (net), and the name of the columns of net that we want to use.

# Run wmean
sample_acts <- run_wmean(mat=counts, net=net, .source='source', .target='target',
                  .mor='weight', times = 100, minsize = 5)
sample_acts
#> # A tibble: 252 × 5
#>    statistic  source   condition            score p_value
#>    <chr>      <chr>    <chr>                <dbl>   <dbl>
#>  1 corr_wmean Androgen PANC1.FOXA2KO.Rep1  4.24     0.26 
#>  2 corr_wmean Androgen PANC1.FOXA2KO.Rep2  6.07     0.140
#>  3 corr_wmean Androgen PANC1.FOXA2KO.Rep3  2.87     0.4  
#>  4 corr_wmean Androgen PANC1.WT.Rep1       7.10     0.100
#>  5 corr_wmean Androgen PANC1.WT.Rep2       3.93     0.28 
#>  6 corr_wmean Androgen PANC1.WT.Rep3       6.59     0.120
#>  7 corr_wmean EGFR     PANC1.FOXA2KO.Rep1  0.0188   0.82 
#>  8 corr_wmean EGFR     PANC1.FOXA2KO.Rep2  0.0224   0.8  
#>  9 corr_wmean EGFR     PANC1.FOXA2KO.Rep3  0.0372   0.66 
#> 10 corr_wmean EGFR     PANC1.WT.Rep1      -0.0520   0.28 
#> # … with 242 more rows

5 Visualization

From the obtained results, we will select the norm_wmean activities and we will observe the obtained activities per sample in a heat-map:

# Transform to wide matrix
sample_acts_mat <- sample_acts %>%
  filter(statistic == 'norm_wmean') %>%
  pivot_wider(id_cols = 'condition', names_from = 'source',
              values_from = 'score') %>%
  column_to_rownames('condition') %>%
  as.matrix()

# Scale per sample
sample_acts_mat <- scale(sample_acts_mat)

# Choose color palette
palette_length = 100
my_color = colorRampPalette(c("Darkblue", "white","red"))(palette_length)

my_breaks <- c(seq(-3, 0, length.out=ceiling(palette_length/2) + 1),
               seq(0.05, 3, length.out=floor(palette_length/2)))

# Plot
pheatmap(sample_acts_mat, border_color = NA, color=my_color, breaks = my_breaks) 

We can observe that WT samples have higher activities for p53 and TGFb than KO. On the other hand, KO show higher activities for MAPK and PI3K.

We can also infer pathway activities from the t-values of the DEGs between KO and WT:

# Run wmean
contrast_acts <- run_wmean(mat=deg, net=net, .source='source', .target='target',
                  .mor='weight', times = 100, minsize = 5)
contrast_acts
#> # A tibble: 42 × 5
#>    statistic  source   condition    score p_value
#>    <chr>      <chr>    <chr>        <dbl>   <dbl>
#>  1 corr_wmean Androgen t         -0.00231  0.0200
#>  2 corr_wmean EGFR     t          0.533    0.1   
#>  3 corr_wmean Estrogen t          0.0550   0.48  
#>  4 corr_wmean Hypoxia  t         -0.00295  0.94  
#>  5 corr_wmean JAK-STAT t          1.90     0.02  
#>  6 corr_wmean MAPK     t          1.99     0.02  
#>  7 corr_wmean NFkB     t          2.04     0.02  
#>  8 corr_wmean PI3K     t          0.962    0.02  
#>  9 corr_wmean TGFb     t         -0.217    0.2   
#> 10 corr_wmean TNFa     t          1.58     0.02  
#> # … with 32 more rows

We select the norm_wmean activities and then we show the changes in activity between KO and WT:

# Filter norm_wmean
f_contrast_acts <- contrast_acts %>%
  filter(statistic == 'norm_wmean')

# Plot
ggplot(f_contrast_acts, aes(x = reorder(source, score), y = score)) + 
    geom_bar(aes(fill = score), stat = "identity") +
    scale_fill_gradient2(low = "darkblue", high = "indianred", 
        mid = "whitesmoke", midpoint = 0) + 
    theme_minimal() +
    theme(axis.title = element_text(face = "bold", size = 12),
        axis.text.x = 
            element_text(angle = 45, hjust = 1, size =10, face= "bold"),
        axis.text.y = element_text(size =10, face= "bold"),
        panel.grid.major = element_blank(), 
        panel.grid.minor = element_blank()) +
    xlab("Pathways")

As observed before, the pathways p53 and TGFb are deactivated in KO when compared to WT, while JAK-STAT and others seem to be activated.

We can further visualize the most responsive genes in each pathway along their t-values to interpret the results. For example, let’s see the genes that are belong to the MAPK pathway:

pathway <- 'MAPK'

df <- net %>%
  filter(source == pathway) %>%
  arrange(target) %>%
  mutate(ID = target, color = "3") %>%
  column_to_rownames('target')
inter <- sort(intersect(rownames(deg),rownames(df)))
df <- df[inter, ]
df['t_value'] <- deg[inter, ]
df <- df %>%
  mutate(color = if_else(weight > 0 & t_value > 0, '1', color)) %>%
  mutate(color = if_else(weight > 0 & t_value < 0, '2', color)) %>%
  mutate(color = if_else(weight < 0 & t_value > 0, '2', color)) %>%
  mutate(color = if_else(weight < 0 & t_value < 0, '1', color))

ggplot(df, aes(x = weight, y = t_value, color = color)) + geom_point() +
  scale_colour_manual(values = c("red","royalblue3","grey")) +
  geom_label_repel(aes(label = ID)) + 
  theme_minimal() +
  theme(legend.position = "none") +
  geom_vline(xintercept = 0, linetype = 'dotted') +
  geom_hline(yintercept = 0, linetype = 'dotted') +
  ggtitle(pathway)
#> Warning: ggrepel: 66 unlabeled data points (too many overlaps). Consider increasing max.overlaps

The pathway seems to be active since the majority of target genes with positive weights have positive t-values (1st quadrant), and the majority of the ones with negative weights have negative t-values (3d quadrant).

6 Session information

#> ─ Session info ───────────────────────────────────────────────────────────────────────────────────────────────────────
#>  setting  value
#>  version  R version 4.2.0 (2022-04-22)
#>  os       Ubuntu 20.04.4 LTS
#>  system   x86_64, linux-gnu
#>  ui       X11
#>  language (EN)
#>  collate  C
#>  ctype    en_US.UTF-8
#>  tz       America/New_York
#>  date     2022-05-15
#>  pandoc   2.5 @ /usr/bin/ (via rmarkdown)
#> 
#> ─ Packages ───────────────────────────────────────────────────────────────────────────────────────────────────────────
#>  package      * version date (UTC) lib source
#>  assertthat     0.2.1   2019-03-21 [2] CRAN (R 4.2.0)
#>  backports      1.4.1   2021-12-13 [2] CRAN (R 4.2.0)
#>  BiocManager    1.30.17 2022-04-22 [2] CRAN (R 4.2.0)
#>  BiocParallel   1.30.2  2022-05-15 [2] Bioconductor
#>  BiocStyle    * 2.24.0  2022-05-15 [2] Bioconductor
#>  bit            4.0.4   2020-08-04 [2] CRAN (R 4.2.0)
#>  bit64          4.0.5   2020-08-30 [2] CRAN (R 4.2.0)
#>  bookdown       0.26    2022-04-15 [2] CRAN (R 4.2.0)
#>  bslib          0.3.1   2021-10-06 [2] CRAN (R 4.2.0)
#>  cellranger     1.1.0   2016-07-27 [2] CRAN (R 4.2.0)
#>  checkmate      2.1.0   2022-04-21 [2] CRAN (R 4.2.0)
#>  cli            3.3.0   2022-04-25 [2] CRAN (R 4.2.0)
#>  colorspace     2.0-3   2022-02-21 [2] CRAN (R 4.2.0)
#>  crayon         1.5.1   2022-03-26 [2] CRAN (R 4.2.0)
#>  curl           4.3.2   2021-06-23 [2] CRAN (R 4.2.0)
#>  data.table     1.14.2  2021-09-27 [2] CRAN (R 4.2.0)
#>  DBI            1.1.2   2021-12-20 [2] CRAN (R 4.2.0)
#>  decoupleR    * 2.2.2   2022-05-15 [1] Bioconductor
#>  digest         0.6.29  2021-12-01 [2] CRAN (R 4.2.0)
#>  dplyr        * 1.0.9   2022-04-28 [2] CRAN (R 4.2.0)
#>  ellipsis       0.3.2   2021-04-29 [2] CRAN (R 4.2.0)
#>  evaluate       0.15    2022-02-18 [2] CRAN (R 4.2.0)
#>  fansi          1.0.3   2022-03-24 [2] CRAN (R 4.2.0)
#>  farver         2.1.0   2021-02-28 [2] CRAN (R 4.2.0)
#>  fastmap        1.1.0   2021-01-25 [2] CRAN (R 4.2.0)
#>  fastmatch      1.1-3   2021-07-23 [2] CRAN (R 4.2.0)
#>  fgsea          1.22.0  2022-05-15 [2] Bioconductor
#>  generics       0.1.2   2022-01-31 [2] CRAN (R 4.2.0)
#>  ggplot2      * 3.3.6   2022-05-03 [2] CRAN (R 4.2.0)
#>  ggrepel      * 0.9.1   2021-01-15 [2] CRAN (R 4.2.0)
#>  glue           1.6.2   2022-02-24 [2] CRAN (R 4.2.0)
#>  gridExtra      2.3     2017-09-09 [2] CRAN (R 4.2.0)
#>  gtable         0.3.0   2019-03-25 [2] CRAN (R 4.2.0)
#>  highr          0.9     2021-04-16 [2] CRAN (R 4.2.0)
#>  hms            1.1.1   2021-09-26 [2] CRAN (R 4.2.0)
#>  htmltools      0.5.2   2021-08-25 [2] CRAN (R 4.2.0)
#>  httr           1.4.3   2022-05-04 [2] CRAN (R 4.2.0)
#>  igraph         1.3.1   2022-04-20 [2] CRAN (R 4.2.0)
#>  jquerylib      0.1.4   2021-04-26 [2] CRAN (R 4.2.0)
#>  jsonlite       1.8.0   2022-02-22 [2] CRAN (R 4.2.0)
#>  knitr          1.39    2022-04-26 [2] CRAN (R 4.2.0)
#>  labeling       0.4.2   2020-10-20 [2] CRAN (R 4.2.0)
#>  later          1.3.0   2021-08-18 [2] CRAN (R 4.2.0)
#>  lattice        0.20-45 2021-09-22 [2] CRAN (R 4.2.0)
#>  lifecycle      1.0.1   2021-09-24 [2] CRAN (R 4.2.0)
#>  logger         0.2.2   2021-10-19 [2] CRAN (R 4.2.0)
#>  lubridate      1.8.0   2021-10-07 [2] CRAN (R 4.2.0)
#>  magick         2.7.3   2021-08-18 [2] CRAN (R 4.2.0)
#>  magrittr       2.0.3   2022-03-30 [2] CRAN (R 4.2.0)
#>  Matrix         1.4-1   2022-03-23 [2] CRAN (R 4.2.0)
#>  munsell        0.5.0   2018-06-12 [2] CRAN (R 4.2.0)
#>  OmnipathR      3.4.0   2022-05-15 [2] Bioconductor
#>  pheatmap     * 1.0.12  2019-01-04 [2] CRAN (R 4.2.0)
#>  pillar         1.7.0   2022-02-01 [2] CRAN (R 4.2.0)
#>  pkgconfig      2.0.3   2019-09-22 [2] CRAN (R 4.2.0)
#>  plyr           1.8.7   2022-03-24 [2] CRAN (R 4.2.0)
#>  prettyunits    1.1.1   2020-01-24 [2] CRAN (R 4.2.0)
#>  progress       1.2.2   2019-05-16 [2] CRAN (R 4.2.0)
#>  purrr          0.3.4   2020-04-17 [2] CRAN (R 4.2.0)
#>  R6             2.5.1   2021-08-19 [2] CRAN (R 4.2.0)
#>  rappdirs       0.3.3   2021-01-31 [2] CRAN (R 4.2.0)
#>  RColorBrewer   1.1-3   2022-04-03 [2] CRAN (R 4.2.0)
#>  Rcpp           1.0.8.3 2022-03-17 [2] CRAN (R 4.2.0)
#>  readr          2.1.2   2022-01-30 [2] CRAN (R 4.2.0)
#>  readxl         1.4.0   2022-03-28 [2] CRAN (R 4.2.0)
#>  RefManageR   * 1.3.0   2020-11-13 [2] CRAN (R 4.2.0)
#>  rlang          1.0.2   2022-03-04 [2] CRAN (R 4.2.0)
#>  rmarkdown      2.14    2022-04-25 [2] CRAN (R 4.2.0)
#>  sass           0.4.1   2022-03-23 [2] CRAN (R 4.2.0)
#>  scales         1.2.0   2022-04-13 [2] CRAN (R 4.2.0)
#>  sessioninfo    1.2.2   2021-12-06 [2] CRAN (R 4.2.0)
#>  stringi        1.7.6   2021-11-29 [2] CRAN (R 4.2.0)
#>  stringr        1.4.0   2019-02-10 [2] CRAN (R 4.2.0)
#>  tibble       * 3.1.7   2022-05-03 [2] CRAN (R 4.2.0)
#>  tidyr        * 1.2.0   2022-02-01 [2] CRAN (R 4.2.0)
#>  tidyselect     1.1.2   2022-02-21 [2] CRAN (R 4.2.0)
#>  tzdb           0.3.0   2022-03-28 [2] CRAN (R 4.2.0)
#>  utf8           1.2.2   2021-07-24 [2] CRAN (R 4.2.0)
#>  vctrs          0.4.1   2022-04-13 [2] CRAN (R 4.2.0)
#>  vroom          1.5.7   2021-11-30 [2] CRAN (R 4.2.0)
#>  withr          2.5.0   2022-03-03 [2] CRAN (R 4.2.0)
#>  xfun           0.31    2022-05-10 [2] CRAN (R 4.2.0)
#>  xml2           1.3.3   2021-11-30 [2] CRAN (R 4.2.0)
#>  yaml           2.3.5   2022-02-21 [2] CRAN (R 4.2.0)
#> 
#>  [1] /tmp/RtmpD6m20a/Rinstb65fd7f3bbb8a
#>  [2] /home/biocbuild/bbs-3.15-bioc/R/library
#> 
#> ──────────────────────────────────────────────────────────────────────────────────────────────────────────────────────