Contents

Compiled date: 2022-04-26

Last edited: 2021-12-12

License: GPL-3

1 Installation

Run the following code to install the Bioconductor version of package.

# install.packages("BiocManager")
BiocManager::install("POMA")

2 Load POMA

library(POMA)

You can also load some additional packages that will be very useful in this vignette.

library(ggplot2)
library(ggraph)
library(plotly)

3 The POMA Workflow

POMA functions can be divided in three sequential well separated blocks: Data Preparation, Pre-processing and Statistical Analysis.

3.1 Data Preparation

The SummarizedExperiment Bioconductor package provides a well defined computational data structures to represent omics experiment data types (Morgan et al. 2020). Since data structures can mean a marked improvement in data analysis, POMA functions use SummarizedExperiment objects from SummarizedExperiment package, allowing the reusability of existing methods for this class and contributing to the improvement of robust and reproducible workflows.

The first step of workflow will be load or create a SummarizedExperiment object. Often, you will have your data stored in separated matrix and/or data frames and you will want to create your SummarizedExperiment object. The PomaSummarizedExperiment function makes this step fast and easy building this SummarizedExperiment object for you.

# create an SummarizedExperiment object from two separated data frames
target <- readr::read_csv("your_target.csv")
features <- readr::read_csv("your_features.csv")

data <- PomaSummarizedExperiment(target = target, features = features)

Alternatively, if your data is already stored in a SummarizedExperiment object, you can skip this step and go directly to the pre-processing step. In this vignette we will use the sample data provided in POMA.

# load example data
data("st000336")
st000336
> class: SummarizedExperiment 
> dim: 31 57 
> metadata(0):
> assays(1): ''
> rownames(31): x1_methylhistidine x3_methylhistidine ... pyruvate
>   succinate
> rowData names(0):
> colnames(57): DMD004.1.U02 DMD005.1.U02 ... DMD167.5.U02 DMD173.1.U02
> colData names(2): group steroids

3.1.1 Brief Description of Example Data

This example data is composed of 57 samples, 31 metabolites, 1 covariate and 2 experimental groups (Controls and DMD) from a targeted LC/MS study.

Duchenne Muscular Dystrophy (DMD) is an X-linked recessive form of muscular dystrophy that affects males via a mutation in the gene for the muscle protein, dystrophin. Progression of the disease results in severe muscle loss, ultimately leading to paralysis and death. Steroid therapy has been a commonly employed method for reducing the severity of symptoms. This study aims to quantify the urine levels of amino acids and organic acids in patients with DMD both with and without steroid treatment. Track the progression of DMD in patients who have provided multiple urine samples.

This data was collected from here.

3.2 Pre Processing

This is a critical point in the workflow because all final statistical results will depend on the decisions made here. Again, this block can be divided in 3 steps: Missing Value Imputation, Normalization and Outlier Detection.

3.2.1 Missing Value Imputation

Often, due to biological and technical reasons, some features can not be identified or quantified in some samples in MS (Armitage et al. 2015). POMA offers 7 different imputation methods to deal with this situation. Just run the following line of code to impute your missings!

imputed <- PomaImpute(st000336, ZerosAsNA = TRUE, RemoveNA = TRUE, cutoff = 20, method = "knn")
imputed
> class: SummarizedExperiment 
> dim: 30 57 
> metadata(0):
> assays(1): ''
> rownames(30): x1_methylhistidine x3_methylhistidine ... pyruvate
>   succinate
> rowData names(0):
> colnames(57): DMD004.1.U02 DMD005.1.U02 ... DMD167.5.U02 DMD173.1.U02
> colData names(2): group steroids

3.2.2 Normalization

The next step of this block is the data normalization. Often, some factors can introduce variability in some types of MS data having a critical influence on the final statistical results, making normalization a key step in the workflow (Berg et al. 2006). Again, POMA offers several methods to normalize the data by running just the following line of code:

normalized <- PomaNorm(imputed, method = "log_pareto")
normalized
> class: SummarizedExperiment 
> dim: 30 57 
> metadata(0):
> assays(1): ''
> rownames(30): x1_methylhistidine x3_methylhistidine ... pyruvate
>   succinate
> rowData names(0):
> colnames(57): DMD004.1.U02 DMD005.1.U02 ... DMD167.5.U02 DMD173.1.U02
> colData names(2): group steroids

3.2.2.1 Normalization effect

Sometimes, you will be interested in how the normalization process affect your data?

To answer this question, POMA offers two exploratory functions, PomaBoxplots and PomaDensity, that can help to understand the normalization process.

PomaBoxplots generates boxplots for all samples or features (depending on the group factor) of a SummarizedExperiment object. Here, we can compare objects before and after normalization step.

PomaBoxplots(imputed, group = "samples", jitter = FALSE) +
  ggtitle("Not Normalized") +
  theme(legend.position = "none") # data before normalization

PomaBoxplots(normalized, group = "samples", jitter = FALSE) +
  ggtitle("Normalized") # data after normalization

On the other hand, PomaDensity shows the distribution of all features before and after the normalization process.

PomaDensity(imputed, group = "features") +
  ggtitle("Not Normalized") +
  theme(legend.position = "none") # data before normalization

PomaDensity(normalized, group = "features") +
  ggtitle("Normalized") # data after normalization

3.2.3 Outlier Detection

Finally, the last step of this block is the Outlier Detection. Outlers are defined as observations that are not concordant with those of the vast majority of the remaining data points. These values can have an enormous influence on the resultant statistical analysis, being a dangerous ground for all required assumptions in the most commonly applied parametric tests in mass spectrometry as well as for all also required assumptions in many regression techniques and predictive modeling approaches. POMA allows the analysis of outliers as well as the possibility to remove them from the analysis using different modulable parameters.

Analyze and remove outliers running the following two lines of code.

PomaOutliers(normalized, do = "analyze")$polygon_plot # to explore

pre_processed <- PomaOutliers(normalized, do = "clean") # to remove outliers
pre_processed
> class: SummarizedExperiment 
> dim: 30 50 
> metadata(0):
> assays(1): ''
> rownames(30): x1_methylhistidine x3_methylhistidine ... pyruvate
>   succinate
> rowData names(0):
> colnames(50): DMD004.1.U02 DMD005.1.U02 ... DMD167.5.U02 DMD173.1.U02
> colData names(2): group steroids

3.3 Statistical Analysis

Once the data have been pre-processed, you can start with the statistical analysis step! POMA offers many different statistical methods and possible combinations to compute. However, in this vignette we will comment only some of the most used.

3.3.1 Univariate Analysis

POMA allows you to perform all of the most used univariate statistical methods in MS by using only one function! PomaUnivariate wrap 4 different univariate methods (ttest, ANOVA and ANCOVA, Wilcoxon test and Kruskal-Wallis Rank Sum Test) that you can perform changing only the “method” argument.

3.3.1.1 T-test

PomaUnivariate(pre_processed, method = "ttest")
> # A tibble: 30 × 9
>    feature                FC diff_means  pvalue pvalueAdj mean_Controls mean_DMD
>    <chr>               <dbl>      <dbl>   <dbl>     <dbl>         <dbl>    <dbl>
>  1 x1_methylhistidine -0.401      0.563 9.30e-8   3.10e-7        -0.402   0.161 
>  2 x3_methylhistidine -0.46       0.613 7.82e-3   9.03e-3        -0.420   0.193 
>  3 alanine            -0.346      0.426 6.62e-4   8.63e-4        -0.317   0.109 
>  4 arginine           -0.558      0.179 4.80e-1   4.80e-1        -0.115   0.0641
>  5 asparagine         -0.367      0.475 1.28e-5   2.40e-5        -0.348   0.128 
>  6 aspartic_acid      -0.35       0.343 3.30e-2   3.54e-2        -0.254   0.0889
>  7 glutamic_acid      -0.274      0.587 1.15e-4   1.82e-4        -0.461   0.126 
>  8 glutamine          -0.477      0.438 5.72e-4   7.80e-4        -0.297   0.141 
>  9 glycine            -0.589      0.412 3.98e-3   4.97e-3        -0.259   0.153 
> 10 histidine          -0.544      0.383 4.13e-4   5.90e-4        -0.248   0.135 
> # … with 20 more rows, and 2 more variables: sd_Controls <dbl>, sd_DMD <dbl>

You can also compute a volcano plot using the T-test results. Note that we’re using the non-normalized object to avoid negative values in our data.

PomaVolcano(imputed, pval = "adjusted")

3.3.1.2 Wilcoxon Test

PomaUnivariate(pre_processed, method = "mann")
> # A tibble: 30 × 9
>    feature                FC diff_means  pvalue pvalueAdj mean_Controls mean_DMD
>    <chr>               <dbl>      <dbl>   <dbl>     <dbl>         <dbl>    <dbl>
>  1 x1_methylhistidine -0.401      0.563 6.21e-5  0.000143        -0.402   0.161 
>  2 x3_methylhistidine -0.46       0.613 9.99e-3  0.0115          -0.420   0.193 
>  3 alanine            -0.346      0.426 2.15e-4  0.000404        -0.317   0.109 
>  4 arginine           -0.558      0.179 1.40e-1  0.145           -0.115   0.0641
>  5 asparagine         -0.367      0.475 2.47e-4  0.000436        -0.348   0.128 
>  6 aspartic_acid      -0.35       0.343 5.71e-3  0.00685         -0.254   0.0889
>  7 glutamic_acid      -0.274      0.587 6.01e-4  0.000859        -0.461   0.126 
>  8 glutamine          -0.477      0.438 2.94e-3  0.00383         -0.297   0.141 
>  9 glycine            -0.589      0.412 4.20e-3  0.00525         -0.259   0.153 
> 10 histidine          -0.544      0.383 2.51e-3  0.00342         -0.248   0.135 
> # … with 20 more rows, and 2 more variables: sd_Controls <dbl>, sd_DMD <dbl>

3.3.2 Limma

Other of the wide used statistical methods in many different omics, such as epigenomics or transcriptomics, is limma (Ritchie et al. 2015). POMA provides an easy use implementation of limma you only have to specify the desired contrast to compute.

PomaLimma(pre_processed, contrast = "Controls-DMD", adjust = "fdr")
> # A tibble: 30 × 7
>    feature             logFC AveExpr     t       P.Value    adj.P.Val     B
>    <chr>               <dbl>   <dbl> <dbl>         <dbl>        <dbl> <dbl>
>  1 tryptophan         -0.775 0.00862 -7.01 0.00000000269 0.0000000807 11.0 
>  2 valine             -0.701 0.0127  -6.63 0.0000000116  0.000000173   9.55
>  3 ornithine          -0.633 0.0337  -6.26 0.0000000479  0.000000479   8.16
>  4 isoleucine         -0.606 0.00438 -5.95 0.000000159   0.00000119    6.99
>  5 lactate            -0.785 0.0184  -5.69 0.000000430   0.00000258    6.02
>  6 pyruvate           -0.624 0.0121  -5.43 0.00000112    0.00000514    5.09
>  7 leucine            -0.616 0.0153  -5.41 0.00000120    0.00000514    5.02
>  8 methionine         -0.552 0.0197  -4.91 0.00000764    0.0000287     3.22
>  9 serine             -0.536 0.0448  -4.76 0.0000132     0.0000440     2.69
> 10 x1_methylhistidine -0.563 0.0373  -4.73 0.0000147     0.0000440     2.59
> # … with 20 more rows

3.3.3 Multivariate Analysis

On the other hand, multivariate analysis implemented in POMA is quite similar to the univariate approaches. PomaMultivariate allows users to compute a PCA, PLS-DA or sPLS-DA by changing only the “method” parameter. This function is based on mixOmics package (Rohart et al. 2017).

3.3.3.1 Principal Component Analysis

poma_pca <- PomaMultivariate(pre_processed, method = "pca")
poma_pca$scoresplot +
  ggtitle("Scores Plot")

3.3.3.2 PLS-DA

poma_plsda <- PomaMultivariate(pre_processed, method = "plsda")
poma_plsda$scoresplot +
  ggtitle("Scores Plot")

poma_plsda$errors_plsda_plot +
  ggtitle("Error Plot")

3.3.4 Correlation Analysis

Often, correlation analysis is used to explore and discover relationships and patterns within our data. PomaCorr provides a flexible and easy way to do that providing a table with all pairwise coorelations in the data, a correlogram and a correlation graph.

poma_cor <- PomaCorr(pre_processed, label_size = 8, coeff = 0.6)
poma_cor$correlations
> # A tibble: 435 × 5
>    feature1      feature2          R   pvalue      FDR
>    <chr>         <chr>         <dbl>    <dbl>    <dbl>
>  1 isoleucine    leucine       0.963 5.03e-29 2.19e-26
>  2 leucine       valine        0.941 3.77e-24 7.11e-22
>  3 fumarate      malate        0.940 4.90e-24 7.11e-22
>  4 isoleucine    valine        0.938 1.07e-23 1.17e-21
>  5 asparagine    threonine     0.907 1.25e-19 1.09e-17
>  6 serine        threonine     0.893 2.83e-18 1.81e-16
>  7 phenylalanine tyrosine      0.893 2.91e-18 1.81e-16
>  8 tryptophan    tyrosine      0.887 9.64e-18 5.24e-16
>  9 leucine       phenylalanine 0.876 8.53e-17 4.13e-15
> 10 asparagine    serine        0.871 2.09e-16 9.07e-15
> # … with 425 more rows
poma_cor$corrplot

poma_cor$graph

Alternatively, if you switch the “corr_type” parameter to “glasso”, this function will compute a Gaussian Graphical Model using the glmnet package (Friedman, Hastie, and Tibshirani 2019).

PomaCorr(pre_processed, corr_type = "glasso", coeff = 0.6)$graph

3.3.5 Lasso, Ridge and Elasticnet

POMA also provides a function to perform a Lasso, Ridge and Elasticnet regression for binary outcomes in a very intuitive and easy way. PomaLasso is based on glmnet package (Friedman, Hastie, and Tibshirani 2010). This function allows you to create a test subset in your data, evaluate the prediction of your models and export the model computed (it could be useful to perform prediction models with MS data). If “ntest” parameter is set to NULL, PomaLasso will use all observations to create the model (useful for feature selection).

# alpha = 1 for Lasso
PomaLasso(pre_processed, alpha = 1, labels = TRUE)$coefficientPlot

3.3.6 Random Forest

Finally, the random forest algorithm is also implemented in POMA. PomaRandForest uses the randomForest package (Liaw and Wiener 2002) to facilitate the implementation of the algorithm and creates automatically both test and train sets to compute and evaluate the resultant models.

poma_rf <- PomaRandForest(pre_processed, ntest = 10, nvar = 10)
poma_rf$error_tree

Resultant random forest model confusion matrix for test set:

poma_rf$confusionMatrix$table
>           Reference
> Prediction 1 2
>          1 2 1
>          2 0 2

Gini index plot for the top 10 predictors:

poma_rf$MeanDecreaseGini_plot

4 Session Information

sessionInfo()
> R version 4.2.0 RC (2022-04-19 r82224)
> Platform: x86_64-pc-linux-gnu (64-bit)
> Running under: Ubuntu 20.04.4 LTS
> 
> Matrix products: default
> BLAS:   /home/biocbuild/bbs-3.15-bioc/R/lib/libRblas.so
> LAPACK: /home/biocbuild/bbs-3.15-bioc/R/lib/libRlapack.so
> 
> 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       
> 
> attached base packages:
> [1] stats     graphics  grDevices utils     datasets  methods   base     
> 
> other attached packages:
> [1] plotly_4.10.0    ggraph_2.0.5     ggplot2_3.3.5    POMA_1.6.0      
> [5] BiocStyle_2.24.0
> 
> loaded via a namespace (and not attached):
>   [1] backports_1.4.1             circlize_0.4.14            
>   [3] plyr_1.8.7                  igraph_1.3.1               
>   [5] lazyeval_0.2.2              splines_4.2.0              
>   [7] gmp_0.6-5                   BiocParallel_1.30.0        
>   [9] listenv_0.8.0               GenomeInfoDb_1.32.0        
>  [11] digest_0.6.29               foreach_1.5.2              
>  [13] htmltools_0.5.2             magick_2.7.3               
>  [15] viridis_0.6.2               fansi_1.0.3                
>  [17] magrittr_2.0.3              cluster_2.1.3              
>  [19] doParallel_1.0.17           limma_3.52.0               
>  [21] graphlayouts_0.8.0          recipes_0.2.0              
>  [23] ComplexHeatmap_2.12.0       globals_0.14.0             
>  [25] gower_1.0.0                 matrixStats_0.62.0         
>  [27] rARPACK_0.11-0              hardhat_0.2.0              
>  [29] colorspace_2.0-3            ggrepel_0.9.1              
>  [31] xfun_0.30                   dplyr_1.0.8                
>  [33] crayon_1.5.1                RCurl_1.98-1.6             
>  [35] jsonlite_1.8.0              impute_1.70.0              
>  [37] survival_3.3-1              iterators_1.0.14           
>  [39] glue_1.6.2                  polyclip_1.10-0            
>  [41] gtable_0.3.0                ipred_0.9-12               
>  [43] zlibbioc_1.42.0             XVector_0.36.0             
>  [45] GetoptLong_1.0.5            DelayedArray_0.22.0        
>  [47] RankProd_3.22.0             future.apply_1.9.0         
>  [49] shape_1.4.6                 Rmpfr_0.8-7                
>  [51] BiocGenerics_0.42.0         scales_1.2.0               
>  [53] DBI_1.1.2                   Rcpp_1.0.8.3               
>  [55] viridisLite_0.4.0           clue_0.3-60                
>  [57] proxy_0.4-26                stats4_4.2.0               
>  [59] lava_1.6.10                 prodlim_2019.11.13         
>  [61] glmnet_4.1-4                httr_1.4.2                 
>  [63] htmlwidgets_1.5.4           RColorBrewer_1.1-3         
>  [65] ellipsis_0.3.2              farver_2.1.0               
>  [67] pkgconfig_2.0.3             nnet_7.3-17                
>  [69] sass_0.4.1                  utf8_1.2.2                 
>  [71] caret_6.0-92                labeling_0.4.2             
>  [73] tidyselect_1.1.2            rlang_1.0.2                
>  [75] reshape2_1.4.4              munsell_0.5.0              
>  [77] tools_4.2.0                 cli_3.3.0                  
>  [79] generics_0.1.2              broom_0.8.0                
>  [81] evaluate_0.15               stringr_1.4.0              
>  [83] fastmap_1.1.0               yaml_2.3.5                 
>  [85] ModelMetrics_1.2.2.2        knitr_1.38                 
>  [87] tidygraph_1.2.1             purrr_0.3.4                
>  [89] randomForest_4.7-1          glasso_1.11                
>  [91] future_1.25.0               nlme_3.1-157               
>  [93] compiler_4.2.0              png_0.1-7                  
>  [95] e1071_1.7-9                 tweenr_1.0.2               
>  [97] tibble_3.1.6                bslib_0.3.1                
>  [99] stringi_1.7.6               highr_0.9                  
> [101] RSpectra_0.16-1             lattice_0.20-45            
> [103] Matrix_1.4-1                vegan_2.6-2                
> [105] permute_0.9-7               vctrs_0.4.1                
> [107] pillar_1.7.0                lifecycle_1.0.1            
> [109] BiocManager_1.30.17         jquerylib_0.1.4            
> [111] GlobalOptions_0.1.2         data.table_1.14.2          
> [113] bitops_1.0-7                corpcor_1.6.10             
> [115] GenomicRanges_1.48.0        R6_2.5.1                   
> [117] bookdown_0.26               gridExtra_2.3              
> [119] IRanges_2.30.0              parallelly_1.31.1          
> [121] codetools_0.2-18            MASS_7.3-57                
> [123] assertthat_0.2.1            SummarizedExperiment_1.26.0
> [125] rjson_0.2.21                withr_2.5.0                
> [127] S4Vectors_0.34.0            GenomeInfoDbData_1.2.8     
> [129] mgcv_1.8-40                 parallel_4.2.0             
> [131] mixOmics_6.20.0             grid_4.2.0                 
> [133] rpart_4.1.16                timeDate_3043.102          
> [135] tidyr_1.2.0                 class_7.3-20               
> [137] rmarkdown_2.14              MatrixGenerics_1.8.0       
> [139] ggforce_0.3.3               pROC_1.18.0                
> [141] Biobase_2.56.0              lubridate_1.8.0            
> [143] ellipse_0.4.2

References

Armitage, Emily Grace, Joanna Godzien, Vanesa Alonso-Herranz, Ángeles López-Gonzálvez, and Coral Barbas. 2015. “Missing Value Imputation Strategies for Metabolomics Data.” Electrophoresis 36 (24): 3050–60.

Berg, Robert A van den, Huub CJ Hoefsloot, Johan A Westerhuis, Age K Smilde, and Mariët J van der Werf. 2006. “Centering, Scaling, and Transformations: Improving the Biological Information Content of Metabolomics Data.” BMC Genomics 7 (1): 142.

Friedman, Jerome, Trevor Hastie, and Rob Tibshirani. 2019. Glasso: Graphical Lasso: Estimation of Gaussian Graphical Models. https://CRAN.R-project.org/package=glasso.

Friedman, Jerome, Trevor Hastie, and Robert Tibshirani. 2010. “Regularization Paths for Generalized Linear Models via Coordinate Descent.” Journal of Statistical Software 33 (1): 1–22. http://www.jstatsoft.org/v33/i01/.

Liaw, Andy, and Matthew Wiener. 2002. “Classification and Regression by randomForest.” R News 2 (3): 18–22. https://CRAN.R-project.org/doc/Rnews/.

Morgan, Martin, Valerie Obenchain, Jim Hester, and Hervé Pagès. 2020. SummarizedExperiment: SummarizedExperiment Container. https://bioconductor.org/packages/SummarizedExperiment.

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. https://doi.org/10.1093/nar/gkv007.

Rohart, Florian, Benoît Gautier, Amrit Singh, and Kim-Anh Lê Cao. 2017. “MixOmics: An R Package for ’Omics Feature Selection and Multiple Data Integration.” PLoS Computational Biology 13 (11): e1005752. http://www.mixOmics.org.