Clustering Mass Spectra from Low Resolution LC-MS/MS Data Using CluMSID

Tobias Depke

November 01, 2022

Introduction

As described in the GC-EI-MS tutorial, CluMSID can also be used to analyse low resolution data – although using low resolution data comes at a cost.

In this example, we will use a similar sample (1uL Pseudomonas aeruginosa PA14 cell extract) as in the General Tutorial, measured with similar chromatography but on a different mass spectrometer, a Bruker amaZon ion trap instrument operated in ESI-(+) mode with auto-MS/MS. In addition to introducing a workflow for low resolution LC-MS/MS data, this example also demonstrates that CluMSID can work with data from different types of mass spectrometers.

Data import

We load the file from the CluMSIDdata package:

library(CluMSID)
library(CluMSIDdata)

lowresfile <- system.file("extdata", 
                        "PA14_amazon_lowres.mzXML",
                        package = "CluMSIDdata")

Data preprocessing

The extraction of spectra works the same way as with high resolution LC-MS/MS data:

ms2list <- extractMS2spectra(lowresfile)
length(ms2list)
#> [1] 1989

Like in the GC-EI-MS example, we have to adjust mz_tolerance to a much higher value compared to high resolution data, while the retention time tolerance can remain unchanged.

featlist <- mergeMS2spectra(ms2list, mz_tolerance = 0.02)
length(featlist)
#> [1] 525

We see that we have similar numbers of spectra as in the General Tutorial, because we tried to keep all parameters except for the mass spectrometer type comparable.

Generation of distance matrix

As we do not have low resolution spectral libraries at hand, we skip the integration of previous knowledge on feature identities in this example and generate a distance matrix right away:

distmat <- distanceMatrix(featlist)

Data exploration

Starting from this distance matrix, we can use all the data exploration functions that CluMSID offers.

When we now make an MDS plot, we learn that the similarity data is very different from the comparable high resolution data:

MDSplot(distmat)
Figure 1: Multidimensional scaling plot as a visualisation of MS2 spectra similarities of the low resolution LC-MS/MS example data set. Black dots signify spectra from unknown metabolites.

Figure 1: Multidimensional scaling plot as a visualisation of MS2 spectra similarities of the low resolution LC-MS/MS example data set. Black dots signify spectra from unknown metabolites.

To get a better overview of the data and the general similarity behaviour, we create a heat map of the distance matrix:

HCplot(distmat, type = "heatmap", 
                cexRow = 0.1, cexCol = 0.1,
                margins = c(6,6))
Figure 2: Symmetric heat map of the distance matrix displaying MS2 spectra similarities of the low resolution LC-MS/MS example data set. along with dendrograms resulting from hierarchical clustering based on the distance matrix. The colour encoding is shown in the top-left insert.

Figure 2: Symmetric heat map of the distance matrix displaying MS2 spectra similarities of the low resolution LC-MS/MS example data set. along with dendrograms resulting from hierarchical clustering based on the distance matrix. The colour encoding is shown in the top-left insert.

We clearly see that the heat map is generally a lot “warmer” than in the General Tutorial (an intuition that is supported by the histogram in the top-left corner), i.e. we have a higher general degree of similarity between spectra. That is not surprising as the m/z information has much fewer levels than in high resolution data and fragments of different sum formula are more likely to have indistinguishable mass-to-charge ratios.

We also see that some more or less compact clusters can be identified. This is easier to inspect in the dendrogram visualisation:

HCplot(distmat, h = 0.8, cex = 0.3)
Figure 3: Circularised dendrogram as a result of agglomerative hierarchical clustering with average linkage as agglomeration criterion based on MS2 spectra similarities of the low resolution LC-MS/MS example data set. Each leaf represents one feature and colours encode cluster affiliation of the features. Leaf labels display feature IDs. Distance from the central point is indicative of the height of the dendrogram.

Figure 3: Circularised dendrogram as a result of agglomerative hierarchical clustering with average linkage as agglomeration criterion based on MS2 spectra similarities of the low resolution LC-MS/MS example data set. Each leaf represents one feature and colours encode cluster affiliation of the features. Leaf labels display feature IDs. Distance from the central point is indicative of the height of the dendrogram.

In conclusion, CluMSID is capable of processing low resolution LC-MS/MS data and if high resolution data is not available, it can be very useful to provide an overview of spectral similarities in low resolution data, thereby helping metabolite annotation if some individual metabolites can be identified by comparison to authentic standards. However, concerning feature annotation, high resolution methods should always be favoured for the many benefits they provide.

Session Info

sessionInfo()
#> R version 4.2.1 (2022-06-23)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 20.04.5 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.16-bioc/R/lib/libRblas.so
#> LAPACK: /home/biocbuild/bbs-3.16-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] stats4    stats     graphics  grDevices utils     datasets  methods  
#> [8] base     
#> 
#> other attached packages:
#>  [1] magrittr_2.0.3      metaMSdata_1.33.0   metaMS_1.34.0      
#>  [4] CAMERA_1.54.0       xcms_3.20.0         MSnbase_2.24.0     
#>  [7] ProtGenerics_1.30.0 S4Vectors_0.36.0    mzR_2.32.0         
#> [10] Rcpp_1.0.9          BiocParallel_1.32.0 Biobase_2.58.0     
#> [13] BiocGenerics_0.44.0 CluMSIDdata_1.13.0  CluMSID_1.14.0     
#> 
#> loaded via a namespace (and not attached):
#>   [1] backports_1.4.1             Hmisc_4.7-1                
#>   [3] plyr_1.8.7                  igraph_1.3.5               
#>   [5] lazyeval_0.2.2              splines_4.2.1              
#>   [7] GenomeInfoDb_1.34.0         ggplot2_3.3.6              
#>   [9] digest_0.6.30               foreach_1.5.2              
#>  [11] htmltools_0.5.3             fansi_1.0.3                
#>  [13] checkmate_2.1.0             cluster_2.1.4              
#>  [15] doParallel_1.0.17           tzdb_0.3.0                 
#>  [17] limma_3.54.0                readr_2.1.3                
#>  [19] sna_2.7                     matrixStats_0.62.0         
#>  [21] vroom_1.6.0                 MsFeatures_1.6.0           
#>  [23] jpeg_0.1-9                  colorspace_2.0-3           
#>  [25] xfun_0.34                   dplyr_1.0.10               
#>  [27] crayon_1.5.2                RCurl_1.98-1.9             
#>  [29] jsonlite_1.8.3              graph_1.76.0               
#>  [31] impute_1.72.0               survival_3.4-0             
#>  [33] iterators_1.0.14            ape_5.6-2                  
#>  [35] glue_1.6.2                  gtable_0.3.1               
#>  [37] zlibbioc_1.44.0             XVector_0.38.0             
#>  [39] DelayedArray_0.24.0         DEoptimR_1.0-11            
#>  [41] scales_1.2.1                vsn_3.66.0                 
#>  [43] DBI_1.1.3                   GGally_2.1.2               
#>  [45] viridisLite_0.4.1           htmlTable_2.4.1            
#>  [47] clue_0.3-62                 archive_1.1.5              
#>  [49] bit_4.0.4                   foreign_0.8-83             
#>  [51] preprocessCore_1.60.0       Formula_1.2-4              
#>  [53] MsCoreUtils_1.10.0          htmlwidgets_1.5.4          
#>  [55] httr_1.4.4                  gplots_3.1.3               
#>  [57] RColorBrewer_1.1-3          ellipsis_0.3.2             
#>  [59] pkgconfig_2.0.3             reshape_0.8.9              
#>  [61] XML_3.99-0.12               farver_2.1.1               
#>  [63] nnet_7.3-18                 sass_0.4.2                 
#>  [65] deldir_1.0-6                utf8_1.2.2                 
#>  [67] labeling_0.4.2              tidyselect_1.2.0           
#>  [69] rlang_1.0.6                 munsell_0.5.0              
#>  [71] tools_4.2.1                 cachem_1.0.6               
#>  [73] cli_3.4.1                   dbscan_1.1-11              
#>  [75] generics_0.1.3              statnet.common_4.7.0       
#>  [77] evaluate_0.17               stringr_1.4.1              
#>  [79] fastmap_1.1.0               mzID_1.36.0                
#>  [81] yaml_2.3.6                  bit64_4.0.5                
#>  [83] knitr_1.40                  robustbase_0.95-0          
#>  [85] caTools_1.18.2              purrr_0.3.5                
#>  [87] RANN_2.6.1                  ncdf4_1.19                 
#>  [89] RBGL_1.74.0                 nlme_3.1-160               
#>  [91] compiler_4.2.1              rstudioapi_0.14            
#>  [93] plotly_4.10.0               png_0.1-7                  
#>  [95] affyio_1.68.0               MassSpecWavelet_1.64.0     
#>  [97] tibble_3.1.8                bslib_0.4.0                
#>  [99] stringi_1.7.8               highr_0.9                  
#> [101] lattice_0.20-45             Matrix_1.5-1               
#> [103] vctrs_0.5.0                 pillar_1.8.1               
#> [105] lifecycle_1.0.3             BiocManager_1.30.19        
#> [107] jquerylib_0.1.4             MALDIquant_1.21            
#> [109] data.table_1.14.4           bitops_1.0-7               
#> [111] GenomicRanges_1.50.0        R6_2.5.1                   
#> [113] latticeExtra_0.6-30         pcaMethods_1.90.0          
#> [115] affy_1.76.0                 network_1.18.0             
#> [117] KernSmooth_2.23-20          gridExtra_2.3              
#> [119] IRanges_2.32.0              codetools_0.2-18           
#> [121] MASS_7.3-58.1               gtools_3.9.3               
#> [123] assertthat_0.2.1            SummarizedExperiment_1.28.0
#> [125] withr_2.5.0                 GenomeInfoDbData_1.2.9     
#> [127] hms_1.1.2                   parallel_4.2.1             
#> [129] grid_4.2.1                  rpart_4.1.19               
#> [131] tidyr_1.2.1                 coda_0.19-4                
#> [133] rmarkdown_2.17              MatrixGenerics_1.10.0      
#> [135] base64enc_0.1-3             interp_1.1-3