--- title: "Creating publication-ready LC-MS data plots" author: "Ossama Edbali" date: "14 December 2025" output: BiocStyle::html_document: toc: true toc_depth: 2 number_sections: true toc_float: true vignette: > %\VignetteIndexEntry{Creating publication-ready LC-MS data plots} %\VignetteEngine{knitr::rmarkdown} %\VignetteEncoding{UTF-8} bibliography: references.bib csl: biomed-central.csl --- ```{r setup, include = FALSE} knitr::opts_knit$set( root.dir = dirname(getwd()) ) knitr::opts_chunk$set( collapse = TRUE, comment = "#>" ) ``` ```{r init, message = FALSE, echo = FALSE, results = "hide"} library(BiocStyle) library(lcmsPlot) ``` # Installation `lcmsPlot` is a Bioconductor package and to install it we have to use [BiocManager](https://cran.r-project.org/web/packages/BiocManager/index.html). ```{r install, eval=FALSE} if (!requireNamespace("BiocManager", quietly = TRUE)) install.packages("BiocManager") BiocManager::install("lcmsPlot") ``` # Load data The dataset used in this vignette originates from a prepared sample spiked with known standards, specifically L-Proline and L-Kynurenine. These compounds are highlighted throughout the vignette and serve as concrete examples to illustrate the different steps of the workflow. To support reproducibility and hands-on exploration, example data files are provided in `inst/extdata/standards-mzml.zip`. Two types of files are included: - An MS1-only file, which is suitable for tasks such as peak detection, feature alignment, and quantitation. - An MS2 file, which additionally enables exploration of fragmentation spectra for compound identification and spectral matching. ```{r load-data, message=FALSE} library(lcmsPlot) exdir <- tempdir() # Unzip the file containing the mzML files utils::unzip( system.file("extdata", "standards-mzml.zip", package = "lcmsPlot"), exdir = exdir ) mzml_dir <- file.path(exdir, "mzml") samples_metadata <- data.frame( sample_id = c( "l-proline-MS1", "l-proline-MS2", "l-kynurenine-MS1", "l-kynurenine-MS2"), sample_name = c(rep("L-Proline", 2), rep("L-Kynurenine", 2)), ms_level = c(1, 2, 1, 2) ) ``` # Overview of the data The following figure illustrates the presence of the two spiked standards, L-Proline and L-Kynurenine, in the MS1-only samples. For each compound, an extracted ion chromatogram is generated using its expected m/z value and retention time, with narrow tolerances to precisely capture the corresponding signal. ```{r data-overview} ms1_samples_metadata <- samples_metadata |> dplyr::filter(ms_level == 1) ms1_raw_files <- file.path( mzml_dir, paste0(ms1_samples_metadata$sample_id, ".mzML")) lcmsPlot(ms1_raw_files, metadata = ms1_samples_metadata) + lp_chromatogram( features = data.frame( sample_id = c("l-proline-MS1", "l-kynurenine-MS1"), mz = c(116.0703793, 209.091824), rt = c(425.74, 365.72) ), ppm = 5, rt_tol = 5 ) + lp_facets(facets = "sample_name", free_x = TRUE, ncol = 2) ``` Next, we plot the chromatograms together with their corresponding MS2 spectra, extracted from the scan closest to the specified retention time for each compound. This allows direct inspection of both chromatographic behavior and fragmentation patterns in a single visualization. ```{r plot-chromatograms-with-spectra} ms2_samples_metadata <- samples_metadata |> dplyr::filter(ms_level == 2) ms2_raw_files <- file.path( mzml_dir, paste0(ms2_samples_metadata$sample_id, ".mzML")) lcmsPlot(ms2_raw_files, metadata = ms2_samples_metadata) + lp_chromatogram( features = data.frame( sample_id = c("l-proline-MS2", "l-kynurenine-MS2"), mz = c(116.0703793, 209.091824), rt = c(425.74, 365.72) ), ppm = 5, rt_tol = 5 ) + lp_spectra( rt = c("l-proline-MS2" = 425.74, "l-kynurenine-MS2" = 365.72), mode = "closest", ms_level = 2) + lp_facets(facets = "sample_name", free_x = TRUE, ncol = 2) ``` # Custom themes You can fully customise the appearance of your plots by modifying the underlying `ggplot2` object. To do this, simply extract the plot using the `lp_get_plot()` function, which returns the raw `ggplot2` object, and then apply any desired theme settings. This gives you complete flexibility to tailor the visualisation to your needs; for example, adjusting fonts, colors, backgrounds, grid lines, or margins. By leveraging custom themes, you can ensure consistency with your project/manuscript's style guidelines, improve readability, or highlight specific aspects of your data. Firstly, define the custom theme using the `ggplot2` `theme` function: ```{r custom-theme} my_custom_theme <- function() { ggplot2::theme( # Title plot.title = ggplot2::element_text( size = 16, margin = ggplot2::margin(b = 15)), # Axis labels axis.title = ggplot2::element_text( size = 11, face = "bold", color = "#1d293d"), # Facet strip strip.text = ggplot2::element_text( color = "#1d293d", size = 9, hjust = 0, margin = ggplot2::margin(t = 3, r = 10, b = 6, l = 10) ), strip.background = ggplot2::element_rect( color = "#90a1b9", fill = "#e2e8f0", linewidth = 0.8), # Box around the panel panel.border = ggplot2::element_rect( color = "#90a1b9", fill = NA, linewidth = 1), # Customize grid lines inside the plot panel.grid.major = ggplot2::element_line( color = "#e2e8f0", linewidth = 0.5), panel.grid.minor = ggplot2::element_line( color = "#e2e8f0", linewidth = 0.5) ) } ``` Generate the plot with `lcmsPlot`, retrieve the underlying plot object using `lp_get_plot()`, and then apply the custom theme: ```{r plot-with-custom-theme} lcmsPlot(ms1_raw_files, metadata = ms1_samples_metadata) + lp_chromatogram( features = data.frame( sample_id = c("l-proline-MS1", "l-kynurenine-MS1"), mz = c(116.0703793, 209.091824), rt = c(425.74, 365.72) ), ppm = 5, rt_tol = 5 ) + lp_facets(facets = "sample_name", free_x = TRUE, ncol = 2) + lp_get_plot() + my_custom_theme() ``` ```{r cleanup, echo=FALSE, results='hide', message=FALSE, warning=FALSE} unlink(mzml_dir, recursive = TRUE, force = TRUE) ``` # Session Info ```{r session-info} sessionInfo() ```