signalomes.Rmd
A key component of the PhosR
package is to construct signalomes. The signalome construction is composed of two main steps: 1) kinase-substrate relationsip scoring and 2) signalome construction. This involves a sequential workflow where the outputs of the first step are used as inputs of the latter step.
In brief, our kinase-substrate relationship scoring method (kinaseSubstrateScore
and kinaseSubstratePred
) prioritises potential kinases that could be responsible for the phosphorylation change of phosphosite on the basis of kinase recognition motif and phosphoproteomic dynamics. Using the kinase-substrate relationships derived from the scoring methods, we reconstruct signalome networks present in the data (Signalomes
) wherin we highlight kinase regulation of discrete modules.
First, we will load the PhosR
package along with few other packages that we will be using in this section of the vignette.
suppressPackageStartupMessages({
library(PhosR)
library(dplyr)
library(ggplot2)
library(GGally)
library(ggpubr)
library(calibrate)
library(network)
})
We will also be needing data containing kinase-substrate annotations from PhosphoSitePlus
, kinase recognition motifs from kinase motifs
, and annotations of kinase families from kinase family
.
As before, we will set up the data by cleaning up the phoshophosite labels and performing RUV normalisation. We will generate the ppe_RUV
matrix as in batch_correction.
data("phospho_L6_ratio_pe")
data("SPSs")
data("PhosphoSitePlus")
##### Run batch correction
ppe <- phospho.L6.ratio.pe
sites = paste(sapply(ppe@GeneSymbol, function(x)x),";",
sapply(ppe@Residue, function(x)x),
sapply(ppe@Site, function(x)x),
";", sep = "")
grps = gsub("_.+", "", colnames(ppe))
design = model.matrix(~ grps - 1)
ctl = which(sites %in% SPSs)
ppe = RUVphospho(ppe, M = design, k = 3, ctl = ctl)
phosphoL6 = ppe@assays@data$normalised
Next, we will filtered for dynamically regulated phosphosites and then standardise the filtered matrix.
# filter for up-regulated phosphosites
phosphoL6.mean <- meanAbundance(phosphoL6, grps = gsub("_.+", "", colnames(phosphoL6)))
aov <- matANOVA(mat=phosphoL6, grps=gsub("_.+", "", colnames(phosphoL6)))
idx <- (aov < 0.05) & (rowSums(phosphoL6.mean > 0.5) > 0)
phosphoL6.reg <- phosphoL6[idx, ,drop = FALSE]
L6.phos.std <- standardise(phosphoL6.reg)
rownames(L6.phos.std) <- paste0(ppe@GeneSymbol, ";", ppe@Residue, ppe@Site, ";")[idx]
We next extract the kinase recognition motifs from each phosphosite.
L6.phos.seq <- ppe@Sequence[idx]
Now that we have all the inputs for kinaseSubstrateScore
and kinaseSubstratePred
ready, we can proceed to the generation of kinase-substrate relationship scores.
L6.matrices <- kinaseSubstrateScore(substrate.list = PhosphoSite.mouse,
mat = L6.phos.std, seqs = L6.phos.seq,
numMotif = 5, numSub = 1, verbose = FALSE)
set.seed(1)
L6.predMat <- kinaseSubstratePred(L6.matrices, top=30, verbose = FALSE)
The signalome construction uses the outputs of kinaseSubstrateScore
and kinaseSubstratePred
functions for the generation of a visualisation of the kinase regulation of discrete regulatory protein modules present in our phosphoproteomic data.
kinaseOI = c("PRKAA1", "AKT1")
Signalomes_results <- Signalomes(KSR=L6.matrices,
predMatrix=L6.predMat,
exprsMat=L6.phos.std,
KOI=kinaseOI)
We can also visualise the relative contribution of each kinase towards the regulation of protein modules by plotting a balloon plot. In the balloon plot, the size of the balloons denote the percentage magnitude of kinase regulation in each module.
### generate palette
my_color_palette <- grDevices::colorRampPalette(RColorBrewer::brewer.pal(8, "Accent"))
kinase_all_color <- my_color_palette(ncol(L6.matrices$combinedScoreMatrix))
names(kinase_all_color) <- colnames(L6.matrices$combinedScoreMatrix)
kinase_signalome_color <- kinase_all_color[colnames(L6.predMat)]
plotSignalomeMap(signalomes = Signalomes_results, color = kinase_signalome_color)
Finally, we can also plot the signalome network that illustrates the connectivity between kinase signalome networks.
plotKinaseNetwork(KSR = L6.matrices, predMatrix = L6.predMat, threshold = 0.9, color = kinase_all_color)
sessionInfo()
#> R version 4.4.1 (2024-06-14)
#> Platform: aarch64-apple-darwin20
#> Running under: macOS Sonoma 14.7
#>
#> Matrix products: default
#> BLAS: /Library/Frameworks/R.framework/Versions/4.4-arm64/Resources/lib/libRblas.0.dylib
#> LAPACK: /Library/Frameworks/R.framework/Versions/4.4-arm64/Resources/lib/libRlapack.dylib; LAPACK version 3.12.0
#>
#> locale:
#> [1] en_US.UTF-8/en_US.UTF-8/en_US.UTF-8/C/en_US.UTF-8/en_US.UTF-8
#>
#> time zone: UTC
#> tzcode source: internal
#>
#> attached base packages:
#> [1] stats graphics grDevices utils datasets methods base
#>
#> other attached packages:
#> [1] network_1.18.2 calibrate_1.7.7 MASS_7.3-60.2 ggpubr_0.6.0
#> [5] GGally_2.2.1 ggplot2_3.5.1 dplyr_1.1.4 PhosR_1.13.1
#>
#> loaded via a namespace (and not attached):
#> [1] bitops_1.0-9 gridExtra_2.3
#> [3] rlang_1.1.4 magrittr_2.0.3
#> [5] matrixStats_1.4.1 e1071_1.7-16
#> [7] compiler_4.4.1 systemfonts_1.1.0
#> [9] vctrs_0.6.5 reshape2_1.4.4
#> [11] stringr_1.5.1 pkgconfig_2.0.3
#> [13] shape_1.4.6.1 crayon_1.5.3
#> [15] fastmap_1.2.0 backports_1.5.0
#> [17] XVector_0.42.0 labeling_0.4.3
#> [19] utf8_1.2.4 rmarkdown_2.28
#> [21] preprocessCore_1.64.0 ragg_1.3.3
#> [23] purrr_1.0.2 xfun_0.48
#> [25] zlibbioc_1.48.2 cachem_1.1.0
#> [27] GenomeInfoDb_1.38.8 jsonlite_1.8.9
#> [29] highr_0.11 DelayedArray_0.28.0
#> [31] broom_1.0.7 R6_2.5.1
#> [33] stringi_1.8.4 bslib_0.8.0
#> [35] RColorBrewer_1.1-3 limma_3.58.1
#> [37] car_3.1-3 GenomicRanges_1.54.1
#> [39] jquerylib_0.1.4 Rcpp_1.0.13
#> [41] SummarizedExperiment_1.32.0 knitr_1.48
#> [43] IRanges_2.36.0 Matrix_1.7-0
#> [45] igraph_2.1.1 tidyselect_1.2.1
#> [47] abind_1.4-8 yaml_2.3.10
#> [49] viridis_0.6.5 lattice_0.22-6
#> [51] tibble_3.2.1 plyr_1.8.9
#> [53] withr_3.0.2 Biobase_2.62.0
#> [55] coda_0.19-4.1 evaluate_1.0.1
#> [57] desc_1.4.3 ggstats_0.7.0
#> [59] proxy_0.4-27 circlize_0.4.16
#> [61] pillar_1.9.0 BiocManager_1.30.25
#> [63] MatrixGenerics_1.14.0 carData_3.0-5
#> [65] stats4_4.4.1 generics_0.1.3
#> [67] RCurl_1.98-1.16 S4Vectors_0.40.2
#> [69] munsell_0.5.1 scales_1.3.0
#> [71] BiocStyle_2.30.0 class_7.3-22
#> [73] glue_1.8.0 pheatmap_1.0.12
#> [75] tools_4.4.1 dendextend_1.18.1
#> [77] ggsignif_0.6.4 fs_1.6.4
#> [79] grid_4.4.1 tidyr_1.3.1
#> [81] colorspace_2.1-1 GenomeInfoDbData_1.2.11
#> [83] Formula_1.2-5 cli_3.6.3
#> [85] ruv_0.9.7.1 textshaping_0.4.0
#> [87] fansi_1.0.6 S4Arrays_1.2.1
#> [89] viridisLite_0.4.2 ggdendro_0.2.0
#> [91] pcaMethods_1.94.0 gtable_0.3.6
#> [93] rstatix_0.7.2 sass_0.4.9
#> [95] digest_0.6.37 BiocGenerics_0.48.1
#> [97] SparseArray_1.2.4 farver_2.1.2
#> [99] htmlwidgets_1.6.4 htmltools_0.5.8.1
#> [101] pkgdown_2.1.1 lifecycle_1.0.4
#> [103] GlobalOptions_0.1.2 statnet.common_4.10.0
#> [105] statmod_1.5.0