Here, we are directly working with the SummarizedExperiment data. For more information on how to create the SummarizedExperiment from a proteomics data set, please refer to the “Get Started” vignette.
The example TMT data set originates from (Biadglegne et al. 2022).
In order to compare the performance of different normalization methods on their ability to detect differentially expressed proteins, we first perform a simple normalization here. For more details about how to normalize data and evaluate the normalization appraoches quantitatively and qualitatively using PRONE, please refer to the “Normalization” vignette.
se_norm <- normalize_se(se, c("IRS_on_RobNorm", "IRS_on_Median", "IRS_on_LoessF", "IRS_on_Quantile"), combination_pattern = "_on_")
#> RobNorm normalization not yet performed. Single RobNorm normalization performed now.
#> RobNorm completed.
#> IRS normalization performed on RobNorm-normalized data completed.
#> Median normalization not yet performed. Single Median normalization performed now.
#> Median completed.
#> IRS normalization performed on Median-normalized data completed.
#> LoessF normalization not yet performed. Single LoessF normalization performed now.
#> LoessF completed.
#> IRS normalization performed on LoessF-normalized data completed.
#> Quantile normalization not yet performed. Single Quantile normalization performed now.
#> Quantile completed.
#> IRS normalization performed on Quantile-normalized data completed.
After having performed normalization and evaluated the different normalization methods via qualitative and quantitative analysis, differential expression analysis can be used to further analyze the differences of the normalization methods.
However before, you need to remove the reference samples in case of a TMT experiment. This can be easily done with the function remove_reference_samples()
.
se_norm <- remove_reference_samples(se_norm)
#> 2 reference samples removed from the SummarizedExperiment object.
First, you need to specify the comparisons you want to perform in DE analysis. For this, the function specifc_comparisons()
was developed which helps to build the right comparison strings. For instance, if you are having a condition column “Condition” in your meta data that is built upon two individual groups, such as “diabetic_3d” referring to diabetic samples at day 3 after operation, you can simply specify sep = "_" and the method only considers comparisons where at least one of the group remains static.
However, you can also just simply create a vector of comparisons to ensure the correct order and handle this vector over to the DE analysis method.
comparisons <- specify_comparisons(se_norm, condition = "Group", sep = NULL, control = NULL)
comparisons <- c("PTB-HC", "TBL-HC", "TBL-PTB", "Rx-PTB")
The function run_DE()
performs the DE analysis on the selected SummarizedExperiment and comparisons. DE analysis can be performed on multiple assays (normalization methods) at once using the already known “ain” parameter. The condition of the SummarizedExperiment object, specified at the beginning, can be used (condition = NULL) or any other column of the meta data can specified. Three methods are available for DE analysis: limma (Ritchie et al. 2015), DEqMS (Zhu et al. 2020), and ROTS (Suomi et al. 2017). The meaning of the other parameters can be extracted from the documentation of the method.
Attention! Due to current issues with DEqMS on Bioconductor and the ongoing submission of PRONE to Bioconductor, the main branch of PRONE no longer includes DEqMS functions. If you need to perform DEqMS, please use the with_DEqMS branch.
de_res <- run_DE(se = se_norm,
comparisons = comparisons,
ain = NULL,
condition = NULL,
DE_method = "limma",
covariate = NULL,
logFC = TRUE,
logFC_up = 1,
logFC_down = -1,
p_adj = TRUE,
alpha = 0.05,
trend = TRUE,
robust = TRUE,
)
#> Condition of SummarizedExperiment used!
#> All assays of the SummarizedExperiment will be used.
#> DE Analysis will not be performed on raw data.
#> log2: DE analysis completed.
#> RobNorm: DE analysis completed.
#> IRS_on_RobNorm: DE analysis completed.
#> Median: DE analysis completed.
#> IRS_on_Median: DE analysis completed.
#> LoessF: DE analysis completed.
#> IRS_on_LoessF: DE analysis completed.
#> Quantile: DE analysis completed.
#> IRS_on_Quantile: DE analysis completed.
When running ROTS you need to specify the R and K parameters. In DEqMS, the DEqMS_PSMs_column need to be specified. This can be any column in the rowData(se)
. For more information on the specific parameters, please refer to the method specific documentation.
If you want to apply other logFC or p-value threshold, there is no need to re-run the DE analysis again. With apply_thresholds()
, you can simply change the threshold values.
To get an overview of the DE results of the different normalization methods, you can visualize the number of significant DE proteins per normalization method in a barplot using plot_overview_DE_bar()
. This plot can be generated in different ways by specifying the “plot_type” parameter:
plot_overview_DE_bar(de_res, ain = NULL, comparisons = comparisons, plot_type = "facet_regulation") + ggplot2::theme(legend.position = "bottom", legend.direction = "vertical")
#> All normalization methods of de_res will be visualized.
You can also just visualize two specific comparisons:
plot_overview_DE_bar(de_res, ain = NULL, comparisons = comparisons[seq_len(2)], plot_type = "facet_comp")
#> All normalization methods of de_res will be visualized.
You can also get an overview of the DE results in form of a heatmap using the plot_overview_DE_tile()
.
plot_overview_DE_tile(de_res)
#> All comparisons of de_res will be visualized.
#> All normalization methods of de_res will be visualized.
Another option is to generate volcano plots for each comparison. The function plot_volcano_DE()
generates a grid of volcano plots for each normalization techniques (facet_norm = TRUE) or for each comparison (facet_comparison = TRUE). A list of volcano plots is returned.
plot_volcano_DE(de_res, ain = NULL, comparisons = comparisons[1], facet_norm = TRUE)
#> All normalization methods of de_res will be visualized.
#> $`PTB-HC`
Furthermore, you can visualize the DE results in form of a heatmap. The function plot_heatmap_DE()
generates a heatmap of the DE results for a specific comparison and normalization method.
plot_heatmap_DE(tuberculosis_TMT_se, tuberculosis_TMT_de_res, ain = NULL, comparison = "PTB-HC", condition = NULL, label_by = NULL, pvalue_column = "adj.P.Val")
#> All assays of the SummarizedExperiment will be used.
#> Label of SummarizedExperiment used!
#> Condition of SummarizedExperiment used!
#> <simpleError in stats::hclust(stats::dist(as.matrix(data))): NA/NaN/Inf in foreign function call (arg 10)>
#> <simpleError in stats::hclust(stats::dist(as.matrix(data))): NA/NaN/Inf in foreign function call (arg 10)>
#> <simpleError in stats::hclust(stats::dist(as.matrix(data))): NA/NaN/Inf in foreign function call (arg 10)>
#> <simpleError in stats::hclust(stats::dist(as.matrix(data))): NA/NaN/Inf in foreign function call (arg 10)>
#> <simpleError in stats::hclust(stats::dist(as.matrix(data))): NA/NaN/Inf in foreign function call (arg 10)>
#> $log2
#>
#> $RobNorm
#>
#> $IRS_on_RobNorm
#>
#> $Median
#>
#> $IRS_on_Median
In case you know some biomarkers of previous publications that should be significantly expressed in your biological context, you can check the number of specified markers that are DE in the different normalization methods.
TODO: plot_coverage_DE_bar
Moreover, you can also intersect the DE results of different normalization methods to see how many DE proteins overlap. You can either plot for each requested comparison an individual upset plot (plot_type = “single”) or stack the number of overlapping DE proteins per comparison (“stacked”). Not only the upset plot(s) are returned, but also a table with the intersections is provided by the plot_upset_DE()
function.
intersections <- plot_upset_DE(de_res, ain = NULL, comparisons = comparisons[seq_len(3)], min_degree = 6, plot_type = "stacked")
#> All normalization methods of de_res will be visualized.
# put legend on top due to very long comparisons
intersections$upset[[2]] <- intersections$upset[[2]] + ggplot2::theme(legend.position = "top", legend.direction = "vertical")
intersections$upset
Additionally, the Jaccard similarity index can be calculated to quantify the similarity of the DE results between the different normalization methods. A individual heatmap can be generated for each comparison (“plot_type = single”), a single heatmap facetted by comparison (“facet_comp”) or a single heatmap taking all comparisons into account (“all”) can be generated.
plot_jaccard_heatmap(de_res, ain = NULL, comparisons = comparisons, plot_type = "all")
#> All normalization methods of de_res will be visualized.
PRONE offers the functionality to extract a consensus set of DEPs based on a selection of normalization methods and a threshold for the number of methods that need to agree on the DE status of a protein. The function get_consensus_DE()
returns a list of consensus DE proteins for either each comparison separately or for all comparisons combined.
utils::sessionInfo()
#> R version 4.4.0 RC (2024-04-16 r86468)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 22.04.4 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.20-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] stats graphics grDevices utils datasets methods base
#>
#> other attached packages:
#> [1] PRONE_0.99.1
#>
#> loaded via a namespace (and not attached):
#> [1] rlang_1.1.4 magrittr_2.0.3
#> [3] GetoptLong_1.0.5 clue_0.3-65
#> [5] matrixStats_1.3.0 compiler_4.4.0
#> [7] png_0.1-8 vctrs_0.6.5
#> [9] reshape2_1.4.4 stringr_1.5.1
#> [11] ProtGenerics_1.37.0 shape_1.4.6.1
#> [13] pkgconfig_2.0.3 crayon_1.5.3
#> [15] fastmap_1.2.0 magick_2.8.3
#> [17] XVector_0.45.0 labeling_0.4.3
#> [19] utf8_1.2.4 rmarkdown_2.27
#> [21] UCSC.utils_1.1.0 preprocessCore_1.67.0
#> [23] purrr_1.0.2 xfun_0.45
#> [25] MultiAssayExperiment_1.31.3 zlibbioc_1.51.1
#> [27] cachem_1.1.0 GenomeInfoDb_1.41.1
#> [29] jsonlite_1.8.8 highr_0.11
#> [31] DelayedArray_0.31.3 BiocParallel_1.39.0
#> [33] parallel_4.4.0 cluster_2.1.6
#> [35] R6_2.5.1 RColorBrewer_1.1-3
#> [37] bslib_0.7.0 stringi_1.8.4
#> [39] ComplexUpset_1.3.3 limma_3.61.2
#> [41] GenomicRanges_1.57.1 jquerylib_0.1.4
#> [43] Rcpp_1.0.12 SummarizedExperiment_1.35.1
#> [45] iterators_1.0.14 knitr_1.47
#> [47] IRanges_2.39.0 splines_4.4.0
#> [49] Matrix_1.7-0 igraph_2.0.3
#> [51] tidyselect_1.2.1 abind_1.4-5
#> [53] yaml_2.3.8 doParallel_1.0.17
#> [55] ggtext_0.1.2 codetools_0.2-20
#> [57] affy_1.83.0 lattice_0.22-6
#> [59] tibble_3.2.1 plyr_1.8.9
#> [61] withr_3.0.0 Biobase_2.65.0
#> [63] evaluate_0.24.0 xml2_1.3.6
#> [65] circlize_0.4.16 pillar_1.9.0
#> [67] affyio_1.75.0 BiocManager_1.30.23
#> [69] MatrixGenerics_1.17.0 DT_0.33
#> [71] foreach_1.5.2 stats4_4.4.0
#> [73] MSnbase_2.31.1 MALDIquant_1.22.2
#> [75] ncdf4_1.22 generics_0.1.3
#> [77] S4Vectors_0.43.0 ggplot2_3.5.1
#> [79] munsell_0.5.1 scales_1.3.0
#> [81] glue_1.7.0 lazyeval_0.2.2
#> [83] tools_4.4.0 data.table_1.15.4
#> [85] mzID_1.43.0 QFeatures_1.15.1
#> [87] vsn_3.73.0 mzR_2.39.0
#> [89] XML_3.99-0.17 Cairo_1.6-2
#> [91] grid_4.4.0 impute_1.79.0
#> [93] tidyr_1.3.1 crosstalk_1.2.1
#> [95] MsCoreUtils_1.17.0 colorspace_2.1-0
#> [97] patchwork_1.2.0 GenomeInfoDbData_1.2.12
#> [99] PSMatch_1.9.0 cli_3.6.3
#> [101] fansi_1.0.6 S4Arrays_1.5.1
#> [103] ComplexHeatmap_2.21.0 dplyr_1.1.4
#> [105] AnnotationFilter_1.29.0 pcaMethods_1.97.0
#> [107] gtable_0.3.5 sass_0.4.9
#> [109] digest_0.6.36 BiocGenerics_0.51.0
#> [111] SparseArray_1.5.10 htmlwidgets_1.6.4
#> [113] rjson_0.2.21 farver_2.1.2
#> [115] htmltools_0.5.8.1 lifecycle_1.0.4
#> [117] httr_1.4.7 GlobalOptions_0.1.2
#> [119] statmod_1.5.0 gridtext_0.1.5
#> [121] MASS_7.3-61
Biadglegne, Fantahun, Johannes R. Schmidt, Kathrin M. Engel, Jörg Lehmann, Robert T. Lehmann, Anja Reinert, Brigitte König, Jürgen Schiller, Stefan Kalkhof, and Ulrich Sack. 2022. “Mycobacterium Tuberculosis Affects Protein and Lipid Content of Circulating Exosomes in Infected Patients Depending on Tuberculosis Disease State.” Biomedicines 10 (4): 783. https://doi.org/10.3390/biomedicines10040783.
Ritchie, Matthew E., Belinda Phipson, Di Wu, Yifang Hu, Charity W. Law, Wei Shi, and Gordon K. Smyth. 2015. “Limma Powers Differential Expression Analyses for Rna-Sequencing and Microarray Studies.” Nucleic Acids Research 43 (7): e47–e47. https://doi.org/10.1093/nar/gkv007.
Suomi, Tomi, Fatemeh Seyednasrollah, Maria K. Jaakkola, Thomas Faux, and Laura L. Elo. 2017. “ROTS: An R Package for Reproducibility-Optimized Statistical Testing.” Edited by Timothée Poisot. PLOS Computational Biology 13 (5): e1005562. https://doi.org/10.1371/journal.pcbi.1005562.
Zhu, Yafeng, Lukas M. Orre, Yan Zhou Tran, Georgios Mermelekas, Henrik J. Johansson, Alina Malyutina, Simon Anders, and Janne Lehtiö. 2020. “Deqms: A Method for Accurate Variance Estimation in Differential Protein Expression Analysis.” Molecular & Cellular Proteomics 19 (6): 1047–57. https://doi.org/10.1074/mcp.TIR119.001646.