Package: MsBackendMassbank
Authors: RforMassSpectrometry Package Maintainer [cre], Michael Witting [aut] (https://orcid.org/0000-0002-1462-4426), Johannes Rainer [aut] (https://orcid.org/0000-0002-6977-7147), Michael Stravs [ctb]
Compiled: Fri Mar 1 18:01:15 2024

1 Introduction

The Spectra package provides a central infrastructure for the handling of Mass Spectrometry (MS) data. The package supports interchangeable use of different backends to import MS data from a variety of sources (such as mzML files). The MsBackendMassbank package allows import and handling MS/MS spectrum data from Massbank. This vignette illustrates the usage of the MsBackendMassbank package to include MassBank data into MS data analysis workflow with the Spectra package in R.

2 Installation

The package can be installed with the BiocManager package. To install BiocManager use install.packages("BiocManager") and, after that, BiocManager::install("MsBackendMassbank") to install this package.

3 Importing MS/MS data from MassBank files

MassBank files (as provided by the Massbank github repository) store normally one library spectrum per file, typically centroided and of MS level 2. In our short example below, we load data from a file containing multiple library spectra per file or from files with each a single spectrum provided with this package. Below we first load all required packages and define the paths to the Massbank files.

library(Spectra)
library(MsBackendMassbank)
fls <- dir(system.file("extdata", package = "MsBackendMassbank"),
           full.names = TRUE, pattern = "txt$")
fls
## [1] "/tmp/RtmpZsfzkO/Rinst264c1562b434d9/MsBackendMassbank/extdata/BSU00001.txt"       
## [2] "/tmp/RtmpZsfzkO/Rinst264c1562b434d9/MsBackendMassbank/extdata/MassBankRecords.txt"
## [3] "/tmp/RtmpZsfzkO/Rinst264c1562b434d9/MsBackendMassbank/extdata/RP000501.txt"       
## [4] "/tmp/RtmpZsfzkO/Rinst264c1562b434d9/MsBackendMassbank/extdata/RP000502.txt"       
## [5] "/tmp/RtmpZsfzkO/Rinst264c1562b434d9/MsBackendMassbank/extdata/RP000503.txt"       
## [6] "/tmp/RtmpZsfzkO/Rinst264c1562b434d9/MsBackendMassbank/extdata/RP000511.txt"       
## [7] "/tmp/RtmpZsfzkO/Rinst264c1562b434d9/MsBackendMassbank/extdata/RP000512.txt"       
## [8] "/tmp/RtmpZsfzkO/Rinst264c1562b434d9/MsBackendMassbank/extdata/RP000513.txt"

MS data can be accessed and analyzed through Spectra objects. Below we create a Spectra with the data from these mgf files. To this end we provide the file names and specify to use a MsBackendMassbank() backend as source to enable data import. First we import from a single file with multiple library spectra.

sps <- Spectra(fls[1],
               source = MsBackendMassbank(),
               backend = MsBackendDataFrame())

With that we have now full access to all imported spectra variables that we list below.

spectraVariables(sps)
##  [1] "msLevel"                 "rtime"                  
##  [3] "acquisitionNum"          "scanIndex"              
##  [5] "dataStorage"             "dataOrigin"             
##  [7] "centroided"              "smoothed"               
##  [9] "polarity"                "precScanNum"            
## [11] "precursorMz"             "precursorIntensity"     
## [13] "precursorCharge"         "collisionEnergy"        
## [15] "isolationWindowLowerMz"  "isolationWindowTargetMz"
## [17] "isolationWindowUpperMz"  "acquistionNum"          
## [19] "accession"               "name"                   
## [21] "smiles"                  "exactmass"              
## [23] "formula"                 "inchi"                  
## [25] "cas"                     "inchikey"               
## [27] "adduct"                  "splash"                 
## [29] "title"

The same is possible with multiple files, each containing a library spectrum.

sps <- Spectra(fls[-1],
               source = MsBackendMassbank(),
               backend = MsBackendDataFrame())

spectraVariables(sps)
##  [1] "msLevel"                 "rtime"                  
##  [3] "acquisitionNum"          "scanIndex"              
##  [5] "dataStorage"             "dataOrigin"             
##  [7] "centroided"              "smoothed"               
##  [9] "polarity"                "precScanNum"            
## [11] "precursorMz"             "precursorIntensity"     
## [13] "precursorCharge"         "collisionEnergy"        
## [15] "isolationWindowLowerMz"  "isolationWindowTargetMz"
## [17] "isolationWindowUpperMz"  "acquistionNum"          
## [19] "accession"               "name"                   
## [21] "smiles"                  "exactmass"              
## [23] "formula"                 "inchi"                  
## [25] "cas"                     "inchikey"               
## [27] "adduct"                  "splash"                 
## [29] "title"

By default the complete metadata is read together with the spectra. This can increase loading time. The different metadata blocks can be skipped which reduces import time. This requires to define an additional data.frame indicating what shall be read.

# create data frame to indicate with metadata blocks shall be read.
metaDataBlocks <- data.frame(metadata = c("ac", "ch", "sp", "ms",
                                          "record", "pk", "comment"),
                             read = rep(TRUE, 7))

sps <- Spectra(fls[-1],
               source = MsBackendMassbank(),
               backeend = MsBackendDataFrame(),
               metaBlock = metaDataBlocks)

# all spectraVariables possible in MassBank are read
spectraVariables(sps)
##   [1] "msLevel"                     "rtime"                      
##   [3] "acquisitionNum"              "scanIndex"                  
##   [5] "dataStorage"                 "dataOrigin"                 
##   [7] "centroided"                  "smoothed"                   
##   [9] "polarity"                    "precScanNum"                
##  [11] "precursorMz"                 "precursorIntensity"         
##  [13] "precursorCharge"             "collisionEnergy"            
##  [15] "isolationWindowLowerMz"      "isolationWindowTargetMz"    
##  [17] "isolationWindowUpperMz"      "acquistionNum"              
##  [19] "accession"                   "name"                       
##  [21] "smiles"                      "exactmass"                  
##  [23] "formula"                     "inchi"                      
##  [25] "cas"                         "inchikey"                   
##  [27] "adduct"                      "splash"                     
##  [29] "title"                       "instrument"                 
##  [31] "instrument_type"             "ms_ms_type"                 
##  [33] "ms_cap_voltage"              "ms_col_gas"                 
##  [35] "ms_desolv_gas_flow"          "ms_desolv_temp"             
##  [37] "ms_frag_mode"                "ms_ionization"              
##  [39] "ms_ionization_energy"        "ms_laser"                   
##  [41] "ms_matrix"                   "ms_mass_accuracy"           
##  [43] "ms_mass_range"               "ms_reagent_gas"             
##  [45] "ms_resolution"               "ms_scan_setting"            
##  [47] "ms_source_temp"              "ms_kinetic_energy"          
##  [49] "ms_electron_current"         "ms_reaction_time"           
##  [51] "chrom_carrier_gas"           "chrom_column"               
##  [53] "chrom_column_temp"           "chrom_column_temp_gradient" 
##  [55] "chrom_flow_gradient"         "chrom_flow_rate"            
##  [57] "chrom_inj_temp"              "chrom_inj_temp_gradient"    
##  [59] "chrom_rti_kovats"            "chrom_rti_lee"              
##  [61] "chrom_rti_naps"              "chrom_rti_uoa"              
##  [63] "chrom_rti_uoa_pred"          "chrom_rt"                   
##  [65] "chrom_rt_uoa_pred"           "chrom_solvent"              
##  [67] "chrom_transfer_temp"         "ims_instrument_type"        
##  [69] "ims_drift_gas"               "ims_drift_time"             
##  [71] "ims_ccs"                     "general_conc"               
##  [73] "compound_class"              "link_cayman"                
##  [75] "link_chebi"                  "link_chembl"                
##  [77] "link_chempdb"                "link_chemspider"            
##  [79] "link_comptox"                "link_hmdb"                  
##  [81] "link_kappaview"              "link_kegg"                  
##  [83] "link_knapsack"               "link_lipidbank"             
##  [85] "link_lipidmaps"              "link_nikkaji"               
##  [87] "link_pubchem"                "link_zinc"                  
##  [89] "scientific_name"             "lineage"                    
##  [91] "link"                        "sample"                     
##  [93] "focus_base_peak"             "focus_derivative_form"      
##  [95] "focus_derivative_mass"       "focus_derivative_type"      
##  [97] "focus_ion_type"              "data_processing_comment"    
##  [99] "data_processing_deprofile"   "data_processing_find"       
## [101] "data_processing_reanalyze"   "data_processing_recalibrate"
## [103] "data_processing_whole"       "deprecated"                 
## [105] "date"                        "authors"                    
## [107] "license"                     "copyright"                  
## [109] "publication"                 "project"                    
## [111] "pknum"                       "comment"
# all NA columns can be dropped
spectraVariables(dropNaSpectraVariables(sps))
##  [1] "msLevel"                   "rtime"                    
##  [3] "acquisitionNum"            "scanIndex"                
##  [5] "dataStorage"               "dataOrigin"               
##  [7] "centroided"                "smoothed"                 
##  [9] "polarity"                  "precScanNum"              
## [11] "precursorMz"               "precursorIntensity"       
## [13] "precursorCharge"           "collisionEnergy"          
## [15] "isolationWindowLowerMz"    "isolationWindowTargetMz"  
## [17] "isolationWindowUpperMz"    "acquistionNum"            
## [19] "accession"                 "name"                     
## [21] "smiles"                    "exactmass"                
## [23] "formula"                   "inchi"                    
## [25] "cas"                       "inchikey"                 
## [27] "adduct"                    "splash"                   
## [29] "title"                     "instrument"               
## [31] "instrument_type"           "ms_ms_type"               
## [33] "ms_frag_mode"              "ms_ionization"            
## [35] "chrom_column"              "chrom_flow_gradient"      
## [37] "chrom_flow_rate"           "chrom_rt"                 
## [39] "chrom_solvent"             "compound_class"           
## [41] "link_chebi"                "link_chemspider"          
## [43] "link_comptox"              "link_kegg"                
## [45] "link_pubchem"              "focus_base_peak"          
## [47] "data_processing_reanalyze" "data_processing_whole"    
## [49] "date"                      "authors"                  
## [51] "license"                   "copyright"                
## [53] "pknum"                     "comment"

Besides default spectra variables, such as msLevel, rtime, precursorMz, we also have additional spectra variables such as the title of each spectrum in the mgf file.

sps$rtime
##                       rtime  rtime  rtime  rtime  rtime  rtime 
## 142.14 142.14 142.14 142.14 142.14 142.14 143.94 143.94 143.94
sps$title
##                                                      
## "L-Tryptophan; LC-ESI-QTOF; MS2; CE: 10; R=; [M+H]+" 
##                                                      
## "L-Tryptophan; LC-ESI-QTOF; MS2; CE: 20; R=; [M+H]+" 
##                                                      
## "L-Tryptophan; LC-ESI-QTOF; MS2; CE: 40; R=; [M+H]+" 
##                                                title 
## "L-Tryptophan; LC-ESI-QTOF; MS2; CE: 10; R=; [M+H]+" 
##                                                title 
## "L-Tryptophan; LC-ESI-QTOF; MS2; CE: 20; R=; [M+H]+" 
##                                                title 
## "L-Tryptophan; LC-ESI-QTOF; MS2; CE: 40; R=; [M+H]+" 
##                                                title 
## "L-Tryptophan; LC-ESI-QTOF; MS2; CE: 10; R=; [M-H]-" 
##                                                title 
## "L-Tryptophan; LC-ESI-QTOF; MS2; CE: 20; R=; [M-H]-" 
##                                                title 
## "L-Tryptophan; LC-ESI-QTOF; MS2; CE: 40; R=; [M-H]-"

In addition we can also access the m/z and intensity values of each spectrum.

mz(sps)
## NumericList of length 9
## [[""]] 74.0233 132.0807 144.0805 146.0598 ... 170.0597 188.0699 205.0965
## [[""]] 74.0232 77.0381 86.0027 91.0539 ... 160.0947 170.0596 171.0625 188.07
## [[""]] 53.0019 53.0383 63.0225 65.0381 ... 158.0817 159.0921 160.0755 170.06
## [["mz"]] 74.0233 132.0807 144.0805 146.0598 ... 170.0597 188.0699 205.0965
## [["mz"]] 74.0232 77.0381 86.0027 91.0539 ... 160.0947 170.0596 171.0625 188.07
## [["mz"]] 53.0019 53.0383 63.0225 65.0381 ... 158.0817 159.0921 160.0755 170.06
## [["mz"]] 72.0095 116.0517 117.0554 159.0935 186.0558 203.0826
## [["mz"]] 72.0094 74.0253 116.0511 117.0539 ... 162.0307 186.0548 203.0818
## [["mz"]] 74.0253 116.0508
intensity(sps)
## NumericList of length 9
## [[""]] 646 980 2114 20052 1248 7628 2036 494048 75708
## [[""]] 10186 142 142 750 138 490 126 ... 14266 1478 1600 16504 1446 109762
## [[""]] 324 184 138 3770 500 800 7214 3238 ... 2802 206 898 162 166 814 250 1840
## [["intensity"]] 646 980 2114 20052 1248 7628 2036 494048 75708
## [["intensity"]] 10186 142 142 750 138 490 ... 14266 1478 1600 16504 1446 109762
## [["intensity"]] 324 184 138 3770 500 800 7214 ... 206 898 162 166 814 250 1840
## [["intensity"]] 150 200 32 232 80 12162
## [["intensity"]] 460 3630 658 74 28 190 44 142 50 52 634
## [["intensity"]] 282 424

When importing a large number of mgf files, setting nonStop = TRUE prevents the call to stop whenever problematic mgf files are encountered.

sps <- Spectra(fls, source = MsBackendMassbank(), nonStop = TRUE)

4 Accessing the MassBank MySQL database

An alternative to the import of the MassBank data from individual text files (which can take a considerable amount of time) is to directly access the MS/MS data in the MassBank MySQL database. For demonstration purposes we are using here a tiny subset of the MassBank data which is stored as a SQLite database within this package.

4.1 Pre-requisites

At present it is not possible to directly connect to the main MassBank production MySQL server, thus, to use the MsBackendMassbankSql backend it is required to install the database locally. The MySQL database dump for each MassBank release can be downloaded from here. This dump could be imported to a local MySQL server.

4.2 Direct access to the MassBank database

To use the MsBackendMassbankSql it is required to first connect to a MassBank database. Below we show the R code which could be used for that - but the actual settings (user name, password, database name, or host) will depend on where and how the MassBank database was installed.

library(RMariaDB)
con <- dbConnect(MariaDB(), host = "localhost", user = "massbank",
                 dbname = "MassBank")

To illustrate the general functionality of this backend we use a tiny subset of the MassBank (release 2020.10) which is provided as an small SQLite database within this package. Below we connect to this database.

library(RSQLite)
con <- dbConnect(SQLite(), system.file("sql", "minimassbank.sqlite",
                                       package = "MsBackendMassbank"))

We next initialize the MsBackendMassbankSql backend which supports direct access to the MassBank in a SQL database and create a Spectra object from that.

mb <- Spectra(con, source = MsBackendMassbankSql())
mb
## MSn data (Spectra) with 70 spectra in a MsBackendMassbankSql backend:
##       msLevel precursorMz  polarity
##     <integer>   <numeric> <integer>
## 1           2         506         0
## 2          NA          NA         1
## 3          NA          NA         0
## 4          NA          NA         1
## 5          NA          NA         0
## ...       ...         ...       ...
## 66          2     185.028         0
## 67          2     455.290         1
## 68          2     253.051         0
## 69          2     358.238         1
## 70          2     256.170         1
##  ... 42 more variables/columns.
##  Use  'spectraVariables' to list all of them.

We can now use this Spectra object to access and use the MassBank data for our analysis. Note that the Spectra object itself does not contain any data from MassBank. Any data will be fetched on demand from the database backend.

To get a listing of all available annotations for each spectrum (the so-called spectra variables) we can use the spectraVariables function.

spectraVariables(mb)
##  [1] "msLevel"                 "rtime"                  
##  [3] "acquisitionNum"          "scanIndex"              
##  [5] "dataStorage"             "dataOrigin"             
##  [7] "centroided"              "smoothed"               
##  [9] "polarity"                "precScanNum"            
## [11] "precursorMz"             "precursorIntensity"     
## [13] "precursorCharge"         "collisionEnergy"        
## [15] "isolationWindowLowerMz"  "isolationWindowTargetMz"
## [17] "isolationWindowUpperMz"  "spectrum_id"            
## [19] "spectrum_name"           "date"                   
## [21] "authors"                 "license"                
## [23] "copyright"               "publication"            
## [25] "splash"                  "compound_id"            
## [27] "adduct"                  "ionization"             
## [29] "ionization_voltage"      "fragmentation_mode"     
## [31] "collision_energy_text"   "instrument"             
## [33] "instrument_type"         "formula"                
## [35] "exactmass"               "smiles"                 
## [37] "inchi"                   "inchikey"               
## [39] "cas"                     "pubchem"                
## [41] "synonym"                 "precursor_mz_text"      
## [43] "compound_name"

Through the MsBackendMassbankSql we can thus access spectra information as well as its annotation.

We can access core spectra variables, such as the MS level with the corresponding function msLevel.

head(msLevel(mb))
## [1]  2 NA NA NA NA  2

Spectra variables can also be accessed with $ and the name of the variable. Thus, MS levels can also be accessed with $msLevel:

head(mb$msLevel)
## [1]  2 NA NA NA NA  2

In addition to spectra variables, we can also get the actual peaks (i.e. m/z and intensity values) with the mz and intensity functions:

mz(mb)
## NumericList of length 70
## [[1]] 146.760803 158.863541 174.988785 ... 470.057434 487.989319 585.88446
## [[2]] 22.99 23.07 23.19 38.98 53.15 60.08 ... 391.18 413.22 414.22 429.16 495.2
## [[3]] 108.099566 152.005378 153.01824 ... 341.031964 409.06729 499.103447
## [[4]] 137.0227 138.025605 139.034602 ... 304.271886 316.303106 352.234228
## [[5]] 112.977759 205.053323 208.041128 ... 449.137152 671.1778 673.200866
## [[6]] 78.893929 183.011719 193.957916 ... 328.091949 408.010193 426.021942
## [[7]] 167.03389 168.040758 333.079965 ... 357.073039 365.040698 373.073473
## [[8]] 195.09167 196.095033
## [[9]] 108.096149 152.004656 153.01824 ... 359.997563 409.065314 491.043775
## [[10]] 1278.12 1279.11 1279.18 1306.21 ... 2064.29 2091.32 2092.3 2136.33
## ...
## <60 more elements>

Note that not all spectra from the database were generated using the same instrumentation. Below we list the number of spectra for each type of instrument.

table(mb$instrument_type)
## 
## LC-ESI-ITFT  LC-ESI-QFT   MALDI-TOF 
##           3          50          17

We next subset the data to all spectra from ions generated by electro spray ionization (ESI).

mb <- mb[mb$ionization == "ESI"]
length(mb)
## [1] 50

As a simple example to illustrate the Spectra functionality we next calculate spectra similarity between one spectrum against all other spectra in the database. To this end we use the compareSpectra function with the normalized dot product as similarity function and allowing 20 ppm difference in m/z between matching peaks

library(MsCoreUtils)
sims <- compareSpectra(mb[11], mb[-11], FUN = ndotproduct, ppm = 40)
max(sims)
## [1] 0.7507467

We plot next a mirror plot for the two best matching spectra.

plotSpectraMirror(mb[11], mb[(which.max(sims) + 1)], ppm = 40)

We can also retrieve the compound information for these two best matching spectra. Note that this compounds function works only with the MsBackendMassbankSql backend as it retrieves the corresponding information from the database’s compound annotation table.

mb_match <- mb[c(11, which.max(sims) + 1)]
compounds(mb_match)
## DataFrame with 2 rows and 10 columns
##   compound_id     formula exactmass                 smiles
##     <integer> <character> <numeric>            <character>
## 1          31    C12H10O2   186.068 COC1=C(C=O)C2=CC=CC=..
## 2          45    C12H10O2   186.068 COC1=CC=C(C=O)C2=CC=..
##                    inchi               inchikey         cas     pubchem
##              <character>            <character> <character> <character>
## 1 InChI=1S/C12H10O2/c1.. YIQGLTKAOHRZOL-UHFFF..   1/12/5392   CID:79352
## 2 InChI=1S/C12H10O2/c1.. MVXMNHYVCLMLDD-UHFFF..  15971-29-6   CID:85217
##                                         synonym                   name
##                                 <CharacterList>            <character>
## 1 2-Methoxy-1-naphthal..,2-methoxynaphthalene.. 2-Methoxy-1-naphthal..
## 2 4-Methoxy-1-Naphthal..,4-methoxynaphthalene.. 4-Methoxy-1-Naphthal..

Note that the MsBackendMassbankSql backend does not support parallel processing because the database connection within the backend can not be shared across parallel processes. Any function on a Spectra object that uses a MsBackendMassbankSql will thus (silently) disable any parallel processing, even if the user might have passed one along to the function using the BPPARAM parameter. In general, the backendBpparam function can be used on any Spectra object to test whether its backend supports the provided parallel processing setup (which might be helpful for developers).

5 Session information

sessionInfo()
## R Under development (unstable) (2024-01-16 r85808)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 22.04.4 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.19-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] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
## [1] MsCoreUtils_1.15.4       RSQLite_2.3.5            MsBackendMassbank_1.11.2
## [4] Spectra_1.13.5           ProtGenerics_1.35.3      BiocParallel_1.37.0     
## [7] S4Vectors_0.41.3         BiocGenerics_0.49.1      BiocStyle_2.31.0        
## 
## loaded via a namespace (and not attached):
##  [1] bit_4.0.5              jsonlite_1.8.8         highr_0.10            
##  [4] compiler_4.4.0         BiocManager_1.30.22    Rcpp_1.0.12           
##  [7] blob_1.2.4             magick_2.8.3           parallel_4.4.0        
## [10] cluster_2.1.6          jquerylib_0.1.4        IRanges_2.37.1        
## [13] yaml_2.3.8             fastmap_1.1.1          R6_2.5.1              
## [16] knitr_1.45             MASS_7.3-60.2          bookdown_0.37         
## [19] DBI_1.2.2              bslib_0.6.1            rlang_1.1.3           
## [22] cachem_1.0.8           xfun_0.42              fs_1.6.3              
## [25] sass_0.4.8             bit64_4.0.5            memoise_2.0.1         
## [28] cli_3.6.2              magrittr_2.0.3         digest_0.6.34         
## [31] MetaboCoreUtils_1.11.2 lifecycle_1.0.4        clue_0.3-65           
## [34] vctrs_0.6.5            evaluate_0.23          codetools_0.2-19      
## [37] rmarkdown_2.25         pkgconfig_2.0.3        tools_4.4.0           
## [40] htmltools_0.5.7