1 Introduction

Visualisation of LC–MS data, and, more broadly, X-MS data, is essential for both quality control and effective communication of results. Diagnostic plots can reveal acquisition issues such as retention-time shifts, inconsistencies in internal standards, and signal drift, while also providing publication-ready figures.

Current LC–MS visualisation approaches primarily include:

  • xcms: Base R plotting functions that address common use cases but offer limited flexibility and visual refinement.
  • Workflow-specific tools (e.g. LCMSQA, MetaboAnalystR, CAMERA): Provide visualisation capabilities that are tightly coupled to their respective workflows, limiting reuse and general applicability.
  • Custom ggplot2-based solutions: Highly flexible, but require substantial familiarity with ggplot2 and can become cumbersome for advanced or routine MS-specific visualisations.

Overall, existing solutions tend to be either overly specialised or insufficiently flexible, making them poorly suited for reusable, cross-platform LC–MS visualisation workflows.

1.1 Objectives

lcmsPlot aims to provide a modern, generic solution for LC-MS visualisation, focused on:

  • Flexible plotting: chromatograms, mass traces, spectra, and more.
  • Building-blocks design: ggplot-like interface for customisable, reproducible plots.
  • Interoperability: works seamlessly with xcms outputs or raw formats (e.g., mzML).
  • Efficiency: scales better to large datasets than current solutions.

The goal is for lcmsPlot to become the community standard for LC-MS data visualisation: versatile for exploration, efficient for large studies, and capable of producing publication-ready figures with minimal effort.

2 Installation

lcmsPlot is a Bioconductor package and to install it we have to use BiocManager.

if (!requireNamespace("BiocManager", quietly = TRUE))
    install.packages("BiocManager")
BiocManager::install("lcmsPlot")

3 Load data

For demonstration purposes we will analyse a small subset of the data from [1] in which the metabolic consequences of the knock-out of the fatty acid amide hydrolase (FAAH) gene in mice was investigated (same data as in the xcms vignette).

library(faahKO)
library(MsExperiment)
library(patchwork)
library(lcmsPlot)

BiocParallel::register(BiocParallel::SerialParam())

sample_indices <- c(1, 2, 5, 6, 7, 8, 11, 12)
cdfs <- dir(
    system.file("cdf", package = "faahKO"),
    full.names = TRUE,
    recursive = TRUE)[sample_indices]

sample_names <- sub(
    basename(cdfs),
    pattern = ".CDF",
    replacement = "",
    fixed = TRUE)

pd <- data.frame(
    sample_name = sample_names,
    sample_group = c(rep("KO", 4), rep("WT", 4)),
    stringsAsFactors = FALSE)

faahko <- readMsExperiment(spectraFiles = cdfs, sampleData = pd)
faahko
#> Object of class MsExperiment 
#>  Spectra: MS1 (10224) 
#>  Experiment data: 8 sample(s)
#>  Sample data links:
#>   - spectra: 8 sample(s) to 10224 element(s).

4 Samples summary

Before examining the detected peaks and their chromatograms, we will first gain an overview of the study samples using a series of plots that highlight different aspects of the data.

4.1 The entry point to plot

The lcmsPlot function is the main starting point for generating visualisations of LC-MS data. It acts as a high-level wrapper that prepares the dataset, configures plotting parameters, and returns a plot object that can be further customised.

p = lcmsPlot(faahko, sample_id_column = "sample_name")
print(p)
#> Object of class lcmsPlotClass 
#>  Data object type: MsExperiment 
#>  Metadata: 8 rows, 6 columns 
#>  Sample ID column: sample_name 
#>  NOTE: No data has been requested to plot.

Through the use of the plotting functions described below, the plot object will display the requested visualisations once it is printed.

4.2 Plot base peak chromatogram (BPC)

Base Peak Chromatograms (BPCs) are useful in LC-MS because they simplify complex data and help quickly identify the most intense signals.

In the example below we plot the BPC for the eight samples:

lcmsPlot(faahko, sample_id_column = "sample_name") +
    lp_chromatogram(ppm = 5, rt_tol = 10, aggregation_fun = "max") +
    lp_facets(facets = "sample_name", ncol = 4)

We can also plot multiple chromatograms in the same plot, overlapped. In the example below, base peak chromatograms are computed for each sample, grouped by sample group, and displayed together in one panel.

lcmsPlot(faahko, sample_id_column = "sample_name") +
    lp_chromatogram(ppm = 5, rt_tol = 10, aggregation_fun = "max") +
    lp_arrange(group_by = "sample_group") +
    lp_labels(title = "Base peak chromatograms", legend = "Sample group") +
    lp_legend(position = "bottom")

4.3 Plot total ion chromatogram (TIC)

The Total Ion Chromatogram (TIC) represents the sum of all ion intensities detected at each time point, providing a comprehensive view of all ions entering the detector, not just the most abundant ones.

lcmsPlot(faahko, sample_id_column = "sample_name") +
    lp_chromatogram(ppm = 5, rt_tol = 10, aggregation_fun = "sum") +
    lp_facets(facets = "sample_name", ncol = 4)

4.4 Total ion current

The Total Ion Current (TIC) reflects the overall amount of ion signal entering the detector at any moment, independent of individual compounds. TIC is widely used for quality control and method assessment. Comparable TIC profiles or distributions across samples indicate consistent injection amounts, stable chromatography, and reliable instrument performance. Conversely, large deviations in TIC can signal issues such as sample loss, ion suppression, contamination, or changes in instrument sensitivity.

Below, we generate a TIC plot for the eight samples under consideration:

lcmsPlot(faahko, sample_id_column = "sample_name") +
    lp_total_ion_current(type = "violin") +
    lp_arrange(group_by = "sample_group") +
    lp_labels(title = "Total ion current", legend = "Sample group")

4.5 Intensity maps

Intensity maps provide a compact, two-dimensional overview of LC-MS data by displaying signal intensity as a function of m/z and retention time. They are a fast and intuitive way to explore global structure in the data and to detect systematic patterns that may not be obvious from chromatograms or summary statistics alone.

In the example below, intensity density maps are generated for two samples (ko15 and wt15) to enable a direct visual comparison. The plot is restricted to an m/z range of 200-600 and a retention time window of 4200-4500 seconds. Density-based visualisation is used to smooth the signal and highlight recurring patterns rather than individual noise peaks.

pparam <- BiocParallel::MulticoreParam(workers = 2)
lcmsPlot(faahko, sample_id_column = "sample_name", BPPARAM = pparam) +
    lp_intensity_map(
        sample_ids = c("ko15", "wt15"),
        mz_range = c(200, 600),
        rt_range = c(4200, 4500),
        density = TRUE) +
    lp_facets(facets = "sample_id", ncol = 2)

5 Plot chromatograms

Before plotting chromatograms, we first perform peak detection. This step allows detected chromatographic peaks to be highlighted in subsequent visualisations, while still preserving the ability to plot chromatograms directly from the raw MS data when desired.

Peak detection is carried out using the centWave algorithm, which is well suited for high-resolution LC-MS data.

cwp <- xcms::CentWaveParam(
    peakwidth = c(20, 80),
    noise = 5000,
    prefilter = c(6, 5000))
faahko <- xcms::findChromPeaks(faahko, param = cwp)

Next, we plot extracted ion chromatograms (XICs) for a feature detected by xcms across all samples. The feature is defined by its m/z and retention time limits, and the corresponding chromatographic signal is extracted from each sample.

In this example, we specify the feature manually using its m/z and RT ranges, which correspond to a peak previously detected by xcms. The resulting chromatograms are arranged by sample group to facilitate comparison across experimental conditions.

lcmsPlot(faahko, sample_id_column = "sample_name") +
    lp_chromatogram(features = rbind(c(
        mzmin = 334.9,
        mzmax = 335.1,
        rtmin = 2700,
        rtmax = 2900))) +
    lp_arrange(group_by = "sample_group") +
    lp_labels(legend = "Sample group")

The above plot can be faceted on metadata factors; in this case we facet on feature_id and sample_id:

lcmsPlot(faahko, sample_id_column = "sample_name") +
    lp_chromatogram(
        features = rbind(c(
            mzmin = 334.9,
            mzmax = 335.1,
            rtmin = 2700,
            rtmax = 2900)),
        highlight_peaks = TRUE,
        highlight_peaks_color = "#121212",
        highlight_apices = list(top_n = 1)
    ) +
    lp_facets(facets = c('feature_id', 'sample_id'), ncol = 4) +
    lp_labels(legend = "Sample") +
    lp_legend(position = "bottom")

Moreover, with the highlight_apices option, the apex of each detected chromatographic peak is highlighted in the plot and annotated with its corresponding retention time value. This makes it easier to identify peak maxima and to compare peak positions across samples.

On the other hand, gridded plots (using the lp_grid function) arrange panels in a matrix defined by two variables. This layout enables structured, two-dimensional comparisons, making it easier to examine how chromatographic patterns vary simultaneously across both factors.

lcmsPlot(faahko, sample_id_column = 'sample_name') +
    lp_chromatogram(
        features = rbind(
            c(mzmin = 334.9, mzmax = 335.1, rtmin = 2700, rtmax = 2900),
            c(mzmin = 278.99721, mzmax = 279.00279, rtmin = 2740, rtmax = 2840)
        ),
        highlight_peaks = TRUE,
        highlight_peaks_color = '#f00') +
    lp_grid(rows = 'feature_id', cols = 'sample_id', free_y = TRUE)

In some cases, you may want to highlight specific regions of an extracted ion chromatogram, for example, to emphasize a peak that has not been automatically detected. This can be achieved using the lp_rt_line function, which allows you to add vertical lines at specified retention times, visually marking areas of interest in the chromatogram.

lcmsPlot(faahko, sample_id_column = 'sample_name') +
    lp_chromatogram(features = rbind(c(
        mzmin = 334.9,
        mzmax = 335.1,
        rtmin = 2700,
        rtmax = 2900))) +
    lp_facets(facets = 'sample_id', ncol = 4) +
    lp_rt_line(intercept = 2800, line_type = 'solid', color = 'red') +
    lp_rt_line(intercept = 2750, line_type = 'dashed', color = 'blue') +
    lp_rt_line(intercept = 2850, line_type = 'dashed', color = 'blue')

6 Plot retention time alignment effect

Retention times can vary slightly across LC-MS runs due to instrumental drift or other experimental factors. To accurately compare chromatographic features across samples, it is important to align retention times so that the same feature appears at the same RT in all samples.

Here, we perform retention time alignment using the Obiwarp algorithm:

faahko <- xcms::adjustRtime(faahko, param = xcms::ObiwarpParam(binSize = 0.6))

Now, we can visualise the BPC for the raw and RT-corrected datasets:

lcmsPlot(faahko, sample_id_column = "sample_name") +
    lp_chromatogram(
        ppm = 5,
        rt_tol = 10,
        aggregation_fun = "max",
        rt_type = "both") +
    lp_arrange(group_by = "sample_group") +
    lp_facets("rt_adjusted", ncol = 1) +
    lp_labels(legend = "Sample group") +
    lp_legend(position = "bottom")

The following plot visualises the retention time correction applied during alignment: the x-axis shows the adjusted retention time (rt_adj), while the y-axis displays the difference between the raw and adjusted retention times (rt_raw - rt_adj). This highlights how much each feature’s retention time was shifted to achieve alignment across samples. Smooth, gradual trends in the plot indicate consistent chromatographic drift correction, while abrupt changes may suggest anomalies. This diagnostic view helps you assess the quality and extent of retention time alignment in LC-MS data processing.

lcmsPlot(faahko, sample_id_column = "sample_name") +
    lp_rt_diff_plot() +
    lp_arrange(group_by = "sample_group") +
    lp_labels(legend = "Sample group")

7 Plot grouped peaks across samples (features)

The next step in LC–MS preprocessing is grouping chromatographic peaks across samples. This step groups peaks that likely originate from the same compound but were detected in different samples.

Here, peak grouping is performed using the peak density approach. Peaks are grouped based on their proximity in m/z and retention time, while also considering their presence across samples. The sampleGroups argument specifies the experimental grouping, minFraction defines the minimum fraction of samples within a group in which a peak must be present, and bw controls the retention time bandwidth used for grouping.

pdp <- xcms::PeakDensityParam(
    sampleGroups = sampleData(faahko)$sample_group,
    minFraction = 0.4,
    bw = 30)
faahko <- xcms::groupChromPeaks(faahko, param = pdp)

Each resulting group of peaks is assigned a unique feature ID, defined by a consensus m/z value and retention time. These feature IDs can then be used to automatically extract and plot chromatograms for individual features across all samples, enabling direct visual comparison of aligned and grouped peaks:

lcmsPlot(faahko, sample_id_column = 'sample_name') +
    lp_chromatogram(
        features = c('M205T2792', 'M207T2719', 'M279T2791', 'M284T2719'),
        ppm = 10,
        rt_tol = 80,
        highlight_peaks = TRUE,
        highlight_peaks_factor = "sample_group") +
    lp_arrange(group_by = 'sample_group') +
    lp_facets(facets = 'feature_id', ncol = 2, free_x = TRUE, free_y = TRUE) +
    lp_labels(title = "Four selected features", legend = "Sample") +
    lp_legend(position = "bottom")

8 Plot chromatograms with mass traces

Until now, we have primarily focused on chromatograms to understand the abundance and retention behavior of features across samples. Chromatograms indicate when a compound elutes and how strongly it is detected, but they do not reveal whether the measured mass values are stable and consistent across the chromatographic peak or between samples.

This is where the mass trace (shown in the bottom panel) becomes valuable. Mass traces display how the measured m/z values evolve over the retention-time profile of a peak. They allow assessment of mass accuracy, mass drift, and signal coherence within a peak, helping to determine whether a feature represents a single, well-defined ion or a mixture of overlapping signals. Consistent and narrow mass traces across samples are indicative of high-quality, reproducible features, whereas variability may point to co-eluting compounds, interference, or measurement instability.

lcmsPlot(faahko, sample_id_column = 'sample_name') +
    lp_chromatogram(
        features = 'M205T2792',
        sample_ids = c('ko15', 'ko16', 'wt15', 'wt16'),
        ppm = 10,
        rt_tol = 50,
        highlight_peaks = TRUE) +
    lp_mass_trace() +
    lp_arrange(group_by = 'sample_id') +
    lp_labels(title = "Feature M205T2792", legend = "Sample") +
    lp_legend(position = "bottom")

As for plain chromatograms, we can facet a combination of chromatograms + mass traces on a factor, such as feature_id:

lcmsPlot(faahko, sample_id_column = 'sample_name') +
    lp_chromatogram(
        features = c('M283T3890', 'M205T2792'),
        sample_ids = c('ko15', 'ko16', 'wt15', 'wt16'),
        ppm = 10,
        rt_tol = 50,
        highlight_peaks = TRUE) +
    lp_mass_trace() +
    lp_arrange(group_by = "sample_id") +
    lp_facets(facets = 'feature_id', ncol = 2) +
    lp_labels(legend = "Sample") +
    lp_legend(position = "bottom")

Or multiple factors (feature_id and sample_id):

lcmsPlot(faahko, sample_id_column = 'sample_name') +
    lp_chromatogram(
        features = c('M283T3890', 'M205T2792'),
        sample_ids = c('ko15', 'ko16', 'wt15', 'wt16'),
        ppm = 10,
        rt_tol = 50,
        highlight_peaks = TRUE,
        highlight_peaks_color = '#343434') +
    lp_mass_trace() +
    lp_facets(facets = c('feature_id', 'sample_id'))

9 Plot chromatograms from raw files

Until now, we have been using an MsExperiment object as the data source for plotting. However, lcmsPlot also supports direct extraction of chromatograms from raw MS files, such as those in the mzML format. This can be useful when you want to quickly inspect raw data without first creating a full experiment object or running preprocessing steps.

In the example below, chromatograms are extracted directly from a subset of raw files. Each file is treated as an individual sample and plotted accordingly.

raw_files <- sampleData(faahko)$spectraOrigin[1:4]

lcmsPlot(raw_files) +
    lp_chromatogram(features = rbind(c(
        mzmin = 334.9,
        mzmax = 335.1,
        rtmin = 2700,
        rtmax = 2900))) +
    lp_arrange(group_by = "sample_id") +
    lp_legend(position = "bottom") +
    lp_labels(legend = "Sample")

By default, file names are used to identify samples. To enable more informative grouping, labeling, and faceting, custom metadata can be attached to the raw files and referenced during plotting:

raw_files <- sampleData(faahko)$spectraOrigin[c(1, 2, 5, 6)]

metadata <- data.frame(
    sample_type = c("KO", "KO", "WT", "WT"),
    sample_id = c("S1", "S2", "S3", "S4")
)

lcmsPlot(raw_files, sample_id_column = "sample_id", metadata = metadata) +
    lp_chromatogram(
        features = rbind(
            c(mzmin = 334.9, mzmax = 335.1, rtmin = 2700, rtmax = 2900),
            c(mzmin = 278.99721, mzmax = 279.00279, rtmin = 2740, rtmax = 2840)
        )
    ) +
    lp_grid(rows = 'feature_id', cols = 'sample_type', free_y = TRUE)

10 Plot spectra

While chromatograms describe how signal intensity changes over retention time, mass spectra provide critical compositional insight by revealing which ions contribute to a feature at a specific point in time. Spectra allow inspection of isotopic patterns, adducts, and potential interferences, and are therefore essential for feature validation and annotation.

To include spectra in the visualisation, we use the lp_spectra function. This function extracts and displays mass spectra based on a specified scan selection mode, which determines which scans are used to generate the spectrum.

10.1 Select closest scan to a specified RT

To extract and plot the mass spectrum from the scan closest to a specified retention time, you can use mode = "closest". This is particularly useful when you want to inspect the spectral composition at or near the apex of a chromatographic peak.

In the example below, a chromatogram is first extracted for a single feature and sample. The detected peaks are highlighted, and the mass spectrum corresponding to the scan closest to the specified RT (3860 seconds) is displayed below the chromatogram.

lcmsPlot(faahko, sample_id_column = 'sample_name') +
    lp_chromatogram(
        features = 'M283T3890',
        sample_ids = 'ko15',
        ppm = 5,
        rt_tol = 50,
        highlight_peaks = TRUE) +
    lp_spectra(rt = 3860, mode = "closest", ms_level = 1) +
    lp_labels(title = "Feature M283T3890", legend = "Sample") +
    lp_legend(position = "bottom")

10.2 Select scan closest to a detected peak apex

To extract and plot the mass spectrum from the scan closest to the apex of a detected chromatographic peak, you can use mode = "closest_apex". This mode automatically identifies the retention time of the peak maximum and selects the nearest corresponding scan, removing the need to manually specify an RT value.

In the example below, chromatograms are extracted for two features in a single sample, with detected peaks highlighted. For each feature, the MS1 spectrum corresponding to the closest scan at the peak apex is displayed, allowing direct inspection of the ions contributing to each peak.

lcmsPlot(faahko, sample_id_column = 'sample_name') +
    lp_chromatogram(
        features = rbind(
            c(mzmin = 362.0982, mzmax = 362.1018, rtmin = 3130, rtmax = 3200),
            c(mzmin = 357.1946, mzmax = 357.2054, rtmin = 3800, rtmax = 3900)
        ),
        sample_ids = 'ko15',
        highlight_peaks = TRUE) +
    lp_spectra(mode = "closest_apex", ms_level = 1) +
    lp_facets(facets = "feature_id", ncol = 2) +
    lp_labels(legend = "Sample") +
    lp_legend(position = "bottom") +
    lp_layout(design = "C\nS\nS")

10.3 Select multiple scans across a detected peak

To inspect how the spectral composition changes across the duration of a detected chromatographic peak, multiple scans can be selected within a retention time interval using mode = "across_peak". This mode extracts spectra at regular intervals spanning the peak, enabling assessment of signal stability, co-eluting ions, or changes in ion dominance from peak start to end.

In the example below, a chromatogram is extracted for a single feature in one sample, with the detected peak highlighted. MS1 spectra are then extracted at regular RT intervals (every 10 seconds) across the peak and displayed sequentially beneath the chromatogram.

lcmsPlot(faahko, sample_id_column = 'sample_name') +
    lp_chromatogram(
        features = "M205T2792",
        sample_ids = 'ko15',
        ppm = 5,
        rt_tol = 50,
        highlight_peaks = TRUE) +
    lp_spectra(mode = "across_peak", interval = 10, ms_level = 1) +
    lp_labels(legend = "Sample") +
    lp_legend(position = "bottom") +
    lp_layout(design = "C\nS\nS\nS\nS\nS\nS\nS")

10.4 Plot standalone spectra

So far, spectra have been displayed alongside chromatograms to provide chromatographic context. However, for detailed inspection of spectral features, such as isotope patterns, adducts, or relative ion intensities, it is often preferable to plot spectra on their own, without the accompanying chromatogram.

In the example below, MS1 spectra are extracted from two samples at the scan closest to a specified retention time and displayed together for direct comparison.

lcmsPlot(faahko, sample_id_column = 'sample_name') +
    lp_spectra(
        sample_ids = c('ko15', 'wt15'),
        rt = 2740,
        mode = "closest",
        ms_level = 1) +
    lp_legend(position = "bottom")

11 Interoperability

lcmsPlot is designed to integrate smoothly with results generated by commonly used LC-MS data processing software. This allows users to visualise processed results directly in R, without the need to manually reformat outputs or reprocess raw data. Below we illustrate interoperability with three widely used tools: Compound Discoverer, MZmine, and MS-DIAL.

Note: These examples are for illustration only and require real data to function properly.

11.1 Compound Discoverer

In the example below, a Compound Discoverer results file is loaded and used to extract chromatograms for a subset of compounds of interest. These chromatograms are arranged by sample and compound name, making it easy to compare signal profiles across the dataset. Highlighting detected peaks helps relate the visualised chromatograms back to the original feature detection performed by Compound Discoverer.

lcmsPlot("example_cd.cdResult") +
    lp_compound_discoverer(
        compounds_query = 'name %in% c("Proline", "Betaine")',
        rt_extend = 5
    ) +
    lp_chromatogram(highlight_peaks = TRUE) +
    lp_grid(rows = "sample_id", cols = "name", free_x = TRUE) +
    lp_labels(title = "Compound Discoverer example", legend = "Sample") +
    lp_legend(position = "bottom")

11.2 MZmine

lcmsPlot also supports feature lists generated by MZmine (version 2 and above). Feature tables, raw data files, and optional sample metadata can be combined into a single dataset object that serves as the basis for visualisation.

mzmine_dataset_example <- MZmineFeatureListsSource(
    feature_lists_paths = c(
        "path/to/feature_list_1.csv",
        "path/to/feature_list_2.csv"
    ),
    sample_paths = c(
        "path/to/sample_path_1.mzML",
        "path/to/sample_path_2.mzML"
    ),
    metadata_path = "path/to/metadata.csv"
)

lcmsPlot(mzmine_dataset_example, sample_id_column = 'filename') +
    lp_chromatogram(
        features = rbind(c(mz = 222.09864, rt = 213.636)),
        sample_ids = c("sample_1", "sample_2"),
        rt_tol = 50,
        ppm = 10,
        highlight_peaks = TRUE,
        highlight_peaks_color = "#894545") +
    lp_facets(facets = "sample_id") +
    lp_labels(title = "MZmine example", legend = "Sample") +
    lp_legend(position = "bottom")

11.3 MS-DIAL

Results generated by MS-DIAL can also be directly integrated into lcmsPlot. Peak tables (corresponding to each sample) exported from MS-DIAL are combined with the corresponding raw data files and optional sample metadata to create a unified dataset for visualisation.

msdial_dataset_example <- MsDialPeaksSource(
    peaks_paths = c(
        "path/to/peaks_1.csv",
        "path/to/peaks_2.csv"
    ),
    sample_paths = c(
        "path/to/sample_path_1.mzML",
        "path/to/sample_path_2.mzML"
    ),
    metadata_path = "path/to/metadata.csv"
)

lcmsPlot(msdial_dataset_example, sample_id_column = 'filename') +
    lp_chromatogram(
        features = rbind(c(mz = 222.09864, rt = 213.636)),
        sample_ids = c("sample_1", "sample_2"),
        rt_tol = 50,
        ppm = 10,
        highlight_peaks = TRUE,
        highlight_peaks_color = "#894545") +
    lp_facets(facets = "sample_id") +
    lp_labels(title = "MS-DIAL example", legend = "Sample") +
    lp_legend(position = "bottom")

12 Session Info

sessionInfo()
#> R Under development (unstable) (2026-01-15 r89304)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.3 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.23-bioc/R/lib/libRblas.so 
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0  LAPACK version 3.12.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] patchwork_1.3.2     MsExperiment_1.13.1 ProtGenerics_1.43.0
#> [4] faahKO_1.51.0       lcmsPlot_0.99.16    BiocStyle_2.39.0   
#> 
#> loaded via a namespace (and not attached):
#>   [1] DBI_1.2.3                   rlang_1.1.7                
#>   [3] magrittr_2.0.4              clue_0.3-66                
#>   [5] MassSpecWavelet_1.77.0      otel_0.2.0                 
#>   [7] matrixStats_1.5.0           compiler_4.6.0             
#>   [9] vctrs_0.7.1                 reshape2_1.4.5             
#>  [11] stringr_1.6.0               crayon_1.5.3               
#>  [13] pkgconfig_2.0.3             MetaboCoreUtils_1.19.2     
#>  [15] fastmap_1.2.0               magick_2.9.0               
#>  [17] XVector_0.51.0              labeling_0.4.3             
#>  [19] rmarkdown_2.30              preprocessCore_1.73.0      
#>  [21] xcms_4.9.1                  tinytex_0.58               
#>  [23] purrr_1.2.1                 xfun_0.56                  
#>  [25] MultiAssayExperiment_1.37.2 cachem_1.1.0               
#>  [27] jsonlite_2.0.0              progress_1.2.3             
#>  [29] DelayedArray_0.37.0         BiocParallel_1.45.0        
#>  [31] prettyunits_1.2.0           parallel_4.6.0             
#>  [33] cluster_2.1.8.2             R6_2.6.1                   
#>  [35] bslib_0.10.0                stringi_1.8.7              
#>  [37] RColorBrewer_1.1-3          limma_3.67.0               
#>  [39] GenomicRanges_1.63.1        jquerylib_0.1.4            
#>  [41] iterators_1.0.14            Rcpp_1.1.1                 
#>  [43] Seqinfo_1.1.0               bookdown_0.46              
#>  [45] SummarizedExperiment_1.41.1 knitr_1.51                 
#>  [47] IRanges_2.45.0              Matrix_1.7-4               
#>  [49] igraph_2.2.2                tidyselect_1.2.1           
#>  [51] dichromat_2.0-0.1           abind_1.4-8                
#>  [53] yaml_2.3.12                 doParallel_1.0.17          
#>  [55] codetools_0.2-20            affy_1.89.0                
#>  [57] lattice_0.22-9              tibble_3.3.1               
#>  [59] plyr_1.8.9                  Biobase_2.71.0             
#>  [61] withr_3.0.2                 S7_0.2.1                   
#>  [63] evaluate_1.0.5              isoband_0.3.0              
#>  [65] Spectra_1.21.2              pillar_1.11.1              
#>  [67] affyio_1.81.0               BiocManager_1.30.27        
#>  [69] MatrixGenerics_1.23.0       foreach_1.5.2              
#>  [71] stats4_4.6.0                MSnbase_2.37.0             
#>  [73] MALDIquant_1.22.3           ncdf4_1.24                 
#>  [75] generics_0.1.4              hms_1.1.4                  
#>  [77] S4Vectors_0.49.0            ggplot2_4.0.2              
#>  [79] scales_1.4.0                glue_1.8.0                 
#>  [81] MsFeatures_1.19.0           lazyeval_0.2.2             
#>  [83] tools_4.6.0                 mzID_1.49.0                
#>  [85] data.table_1.18.2.1         QFeatures_1.21.0           
#>  [87] vsn_3.79.5                  mzR_2.45.0                 
#>  [89] fs_1.6.6                    XML_3.99-0.22              
#>  [91] grid_4.6.0                  impute_1.85.0              
#>  [93] tidyr_1.3.2                 MsCoreUtils_1.23.2         
#>  [95] PSMatch_1.15.1              cli_3.6.5                  
#>  [97] viridisLite_0.4.3           S4Arrays_1.11.1            
#>  [99] dplyr_1.2.0                 AnnotationFilter_1.35.0    
#> [101] pcaMethods_2.3.0            gtable_0.3.6               
#> [103] sass_0.4.10                 digest_0.6.39              
#> [105] BiocGenerics_0.57.0         SparseArray_1.11.10        
#> [107] farver_2.1.2                htmltools_0.5.9            
#> [109] lifecycle_1.0.5             statmod_1.5.1              
#> [111] MASS_7.3-65

1. Saghatelian A, Trauger SA, Want EJ, Hawkins EG, Siuzdak G, Cravatt BF: Assignment of endogenous substrates to enzymes by global metabolite profiling. Biochemistry 2004, 43:14332–9.