Contents

set.seed(123)
library(DepInfeR)
library(tibble)
library(tidyr)
library(ggplot2)
library(BiocParallel)
library(dplyr)

1 Installation

The DepInfeR package is available at bioconductor.org and can be downloaded via Bioc Manager:

if (!requireNamespace("BiocManager", quietly = TRUE)) {
  install.packages("BiocManager")
}
BiocManager::install("DepInfeR")

2 Introduction

DepInfeR is an R package for a computational method that integrates two experimentally accessible input data matrices: the drug sensitivity profiles of cancer cell lines or primary tumors ex-vivo (X), and the drug affinities of a set of proteins (Y), to infer a matrix of molecular protein dependencies of the cancers (ß). DepInfeR deconvolutes the protein inhibition effect on the viability phenotype by using regularized multivariate linear regression. It assigns a “dependence coefficient” to each protein and each sample, and therefore could be used to gain a causal and accurate understanding of functional consequences of genomic aberrations in a heterogeneous disease, as well as to guide the choice of pharmacological intervention for a specific cancer type, sub-type, or an individual patient.

This document provides a walk-through of using the DepInfeR package to infer sample-specific protein dependencies from drug-protein affinity profiling and ex-vivo drug response data.

For more information, please read out preprint on bioRxiv: https://doi.org/10.1101/2022.01.11.475864.

3 Data input

As inputs, DepInfeR requires two types of data: 1) a drug-protein affinity profiling dataset containing affinity values between a set of drugs and a set of proteins, and 2) a drug response dataset containing phenotypic data (e.g. viability measurements) of samples in response to a set of drugs.

In this exemplary walk-through analysis we use 1) drug-protein affinity data from Klaeger et al. 2017, which can be found in the supplementary file of the paper (Table_S1 & Table_S2): https://science.sciencemag.org/content/358/6367/eaan4368/tab-figures-data, and 2) drug response data from the Genomics of Drug Sensitivity in Cancer (GDSC) cancer cell line screening dataset: https://www.cancerrxgene.org/.

A subset of leukemia and breast cancer cell lines was chosen for this analysis. The analyzed cancer types were.

The Her2 status of the breast cancer cell lines was annotated manually.

Genetic background of the cancer cell lines was annotated with their mutational background. The mutational background of the cell lines was retrieved from https://www.cancerrxgene.org/downloads/genetic_features.

This analysis starts with the two data tables that contain the common drugs in both drug-target affinity and drug sensitivity datasets.

data(targetsGDSC, drug_response_GDSC)

A glance at the drug-target affinity table

head(targetsGDSC)
## # A tibble: 6 × 6
##   drugName drugID targetClassification    EC50      Kd targetName
##   <chr>     <dbl> <chr>                  <dbl>   <dbl> <chr>     
## 1 afatinib   1032 High confidence      3493.   3493.   ADK       
## 2 afatinib   1032 High confidence         3.01    2.39 EGFR      
## 3 afatinib   1032 High confidence      3154.   1771.   GAK       
## 4 afatinib   1032 High confidence      3828.   2388.   GRB2      
## 5 afatinib   1032 High confidence      3884.   2811.   MAPK11    
## 6 afatinib   1032 High confidence      2044.   1047.   MAPK14

A glance at the drug sensitivity table

head(drug_response_GDSC)
## # A tibble: 6 × 14
##   Drug      `Drug Id` `Cell line name` `Cosmic sample Id` `TCGA classification`
##   <chr>         <dbl> <chr>                         <dbl> <chr>                
## 1 erlotinib         1 CCRF-CEM                     905952 ALL                  
## 2 erlotinib         1 MOLT-4                       905958 ALL                  
## 3 erlotinib         1 BE-13                        906763 ALL                  
## 4 erlotinib         1 GR-ST                        906877 ALL                  
## 5 erlotinib         1 KARPAS-45                    907272 ALL                  
## 6 erlotinib         1 KE-37                        907277 ALL                  
## # ℹ 9 more variables: Tissue <chr>, `Tissue sub-type` <chr>, IC50 <dbl>,
## #   AUC <dbl>, `Max conc` <dbl>, RMSE <dbl>, `Z score` <dbl>,
## #   `Dataset version` <chr>, synonyms <chr>

4 Pre-processing the drug-protein dataset

We first need to do some pre-processing of both datasets, to turn them into numeric matrices that can be used as inputs for DepInfeR.

Rename BCR to BCR/ABL to avoid confusion with B-cell receptor (BCR)

targetsGDSC <- mutate(targetsGDSC,
  targetName = ifelse(targetName %in% "BCR", "BCR/ABL", targetName)
)

Turn target table into drug-protein affinity matrix

targetMatrix <- filter(
  targetsGDSC,
  targetClassification == "High confidence"
) %>%
  select(drugID, targetName, Kd) %>%
  pivot_wider(names_from = "targetName", values_from = "Kd") %>%
  remove_rownames() %>%
  column_to_rownames("drugID") %>%
  as.matrix()

We provide a function, ProcessTargetResults, for the pre-processing of the drug-protein affinity matrix with Kd values (or optionally other affinity measurement values at roughly normal distribution). This function can perform the following steps:

All steps within this function are optional depending on input data. The transformation steps should be performed if the affinity matrix consists of Kd values. If there are highly correlated features within the affinity matrix, they can be removed using the provided function. The users can also use their own method to preprocess the input drug-target affinity matrix, as long as the distribution of the affinity values is roughly normal and do not have too many highly correlated features.

ProcessTargetResults <- processTarget(targetMatrix,
  KdAsInput = TRUE,
  removeCorrelated = TRUE
)

5 Preparation of the drug response matrix

5.1 Prepare drug response matrix using z-scores

The z-score was chosen as a suitable measurement value for our drug screening response matrix as it acts as normalization for each drug over all cell lines. When working with AUC or IC50 values, a suitable normalization of the values is recommended. In this analysis we used the z-score of the AUC values as the drug response input for DepInfeR.

responseMatrix <- filter(drug_response_GDSC, `Drug Id` %in% targetsGDSC$drugID) %>%
  select(`Drug Id`, `Cell line name`, AUC) %>%
  pivot_wider(names_from = `Cell line name`, values_from = AUC) %>%
  remove_rownames() %>%
  column_to_rownames("Drug Id") %>%
  as.matrix()

5.2 Assessment of missing values

Currently, DepInfeR does not support input data with missing values. Therefore, the missing values in the input datasets need to be properly handled. The entries with missing values can either be removed or imputed.

We firstly check the distribution of our missing values across all cell lines. Below is the plot of remaining cell lines depending on cutoff for the maximum number of missing values per column (cell line).

missTab <- data.frame(
  NA_cutoff = 0,
  remain_celllines = 0, stringsAsFactors = FALSE
)

for (i in 0:138) {
  a <- dim(responseMatrix[, colSums(is.na(responseMatrix)) <= i])[2]
  missTab[i, ] <- c(i, a)
}

ggplot(missTab, aes(x=NA_cutoff, y=remain_celllines)) +
    geom_line() + theme_bw() +
    xlab("Allowed NA values") +
    ylab("Number of remaining cell lines") +
    geom_vline(xintercept = 24, color = "red", linetype = "dashed")
This figure shows the number of cell lines included in the analysis as a function of the percentage of missing values allowed per cell line.

Figure 1: This figure shows the number of cell lines included in the analysis as a function of the percentage of missing values allowed per cell line

In order to keep as many cell lines as possible while reducing the missing data points, we chose one of the elbow points (x=24, shown as the red dashed line) in the plot above as the cut-off for allowed missing values. We will impute the remaining missing values by using the MissForest imputation method.

5.3 Subset for cell lines with less than 24 missing values (based on assessment above)

responseMatrix <- responseMatrix[, colSums(is.na(responseMatrix)) <= 24]

5.4 MissForest imputation

library(missForest)

impRes <- missForest(t(responseMatrix))
imp_missforest <- impRes$ximp

responseMatrix_imputed <- t(imp_missforest)

5.5 Calculate column-wise z-score

responseMatrix_scaled <- t(scale(t(responseMatrix_imputed)))

6 Combine the feature and response matrix for the regression model

In this step, we will also subset the drug-target affinity matrix and the drug response matrix to only keep the drugs present in both matrices.

targetInput <- ProcessTargetResults$targetMatrix

overlappedDrugs <- intersect(
  rownames(responseMatrix_scaled),
  rownames(targetInput)
)

targetInput <- targetInput[overlappedDrugs, ]
responseInput <- responseMatrix_scaled[overlappedDrugs, ]

7 Multivariate model for protein dependence prediction

7.1 Multi-target LASSO model

Perform multivariate LASSO regression based on a drug-protein affinity matrix and a drug response matrix.
In this walk-through, we only use 20 repetitions to save time. In a real application, the number of repetitions can be larger, e.g. repeats = 100, to gain more robust estimations. We will also set a seed for the random number generator to ensure reproducibility. Note that the global random number stream by set.seed() function is ignored by the runLASSORegression() function.

#set up BiocParallel back-end
param <- MulticoreParam(workers = 2, RNGseed = 333)

result <- runLASSORegression(
  TargetMatrix = targetInput,
  ResponseMatrix = responseInput, repeats = 3,
  BPPARAM = param
)

Remove targets that were never selected

useTar <- rowSums(result$coefMat != 0) > 0
result$coefMat <- result$coefMat[useTar, ]

Number of selected targets

nrow(result$coefMat)
## [1] 21

8 Examples of how to interpret and perform downstream analyses on the inferred protein dependence matrix

8.1 Heatmap of protein dependence coefficients

The protein dependence matrix can be visualized in a heatmap. High positive coefficients imply strong reliance of a certain sample on this protein for survival. Proteins with coefficients close to zero are less essential for the survival of cells. Negative coefficients indicate that the viability phenotype benefits from inhibition of the protein.

library(RColorBrewer)
library(pheatmap)

# firstly, we need to load the processed mutation annotation of the cell lines
data("mutation_GDSC")

#setup column annotation and color scheme
annoColor <- list(
  H2O2 = c(`-1` = "red", `0` = "black", `1` = "green"),
  IL.1 = c(`-1` = "red", `0` = "black", `1` = "green"),
  JAK.STAT = c(`-1` = "red", `0` = "black", `1` = "green"),
  MAPK.only = c(`-1` = "red", `0` = "black", `1` = "green"),
  MAPK.PI3K = c(`-1` = "red", `0` = "black"),
  TLR = c(`-1` = "red", `0` = "black", `1` = "green"),
  Wnt = c(`-1` = "red", `0` = "black", `1` = "green"),
  VEGF = c(`-1` = "red", `0` = "black", `1` = "green"),
  PI3K.only = c(`-1` = "red", `0` = "black", `1` = "green"),
  TCGA.classification =
    c(
      ALL = "#BC3C29FF", AML = "#E18727FF",
      DLBC = "#20854EFF", "BRCAHer-" = "#0072B5FF",
      "BRCAHer+" = "#7876B1FF"
    ),
  ARID1A_mut = c(`1` = "black", `0` = "grey80"),
  EP300_mut = c(`1` = "black", `0` = "grey80"),
  PTEN_mut = c(`1` = "black", `0` = "grey80"),
  TP53_mut = c(`1` = "black", `0` = "grey80"),
  PIK3CA_mut = c(`1` = "black", `0` = "grey80"),
  BRCA2_mut = c(`1` = "black", `0` = "grey80"),
  BRCA1_mut = c(`1` = "black", `0` = "grey80"),
  CDH1_mut = c(`1` = "black", `0` = "grey80"),
  FBXW7_mut = c(`1` = "black", `0` = "grey80"),
  NRAS_mut = c(`1` = "black", `0` = "grey80"),
  ASXL1_mut = c(`1` = "black", `0` = "grey80"),
  MLL2_mut = c(`1` = "black", `0` = "grey80"),
  ABL1_trans = c(`1` = "black", `0` = "grey80"),
  missing_value_perc = c(`0` = "white", `25` = "red")
)
importanceTab <- result$coefMat
#change "-" to "." in the column name to match cell line annotation
colnames(importanceTab) <- gsub("-",".", colnames(importanceTab),)
plotTab <- importanceTab
# Row normalization while keeping sign
plotTab_scaled <- scale(t(plotTab), center = FALSE, scale = TRUE)
plotTab <- plotTab_scaled
mutation_GDSC$TCGA.classification[mutation_GDSC$TCGA.classification == "BRCA"] <- "BRCAHer-"
mutation_GDSC$TCGA.classification <-
  factor(mutation_GDSC$TCGA.classification,
    levels = c("ALL", "AML", "DLBC", "BRCAHer-", "BRCAHer+")
  )
pheatmap(plotTab,
  color = colorRampPalette(rev(brewer.pal(n = 7, name = "RdBu")),
    bias = 1.8
  )(100),
  annotation_row = mutation_GDSC,
  annotation_colors = annoColor,
  clustering_method = "ward.D2", scale = "none",
  main = "", fontsize = 9,
  fontsize_row = 7, fontsize_col = 10
)