lcmsPlot is a Bioconductor package and to install it we
have to use BiocManager.
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:
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)
)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.
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.
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)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:
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:
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()sessionInfo()
#> R version 4.5.2 (2025-10-31)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.3 LTS
#>
#> Matrix products: default
#> BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
#> LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
#>
#> locale:
#> [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
#> [3] LC_TIME=en_US.UTF-8 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: Etc/UTC
#> tzcode source: system (glibc)
#>
#> attached base packages:
#> [1] stats graphics grDevices utils datasets methods base
#>
#> other attached packages:
#> [1] lcmsPlot_0.99.16 BiocStyle_2.39.0
#>
#> loaded via a namespace (and not attached):
#> [1] sass_0.4.10 generics_0.1.4 digest_0.6.39
#> [4] magrittr_2.0.4 evaluate_1.0.5 grid_4.5.2
#> [7] RColorBrewer_1.1-3 fastmap_1.2.0 jsonlite_2.0.0
#> [10] ProtGenerics_1.43.0 mzR_2.45.0 BiocManager_1.30.27
#> [13] scales_1.4.0 codetools_0.2-20 jquerylib_0.1.4
#> [16] cli_3.6.5 rlang_1.1.7 Biobase_2.71.0
#> [19] withr_3.0.2 cachem_1.1.0 yaml_2.3.12
#> [22] otel_0.2.0 tools_4.5.2 parallel_4.5.2
#> [25] BiocParallel_1.45.0 dplyr_1.2.0 ncdf4_1.24
#> [28] ggplot2_4.0.2 BiocGenerics_0.57.0 buildtools_1.0.0
#> [31] vctrs_0.7.1 R6_2.6.1 lifecycle_1.0.5
#> [34] pkgconfig_2.0.3 pillar_1.11.1 bslib_0.10.0
#> [37] gtable_0.3.6 glue_1.8.0 Rcpp_1.1.1
#> [40] xfun_0.56 tibble_3.3.1 tidyselect_1.2.1
#> [43] sys_3.4.3 knitr_1.51 farver_2.1.2
#> [46] htmltools_0.5.9 patchwork_1.3.2 rmarkdown_2.30
#> [49] maketools_1.3.2 labeling_0.4.3 compiler_4.5.2
#> [52] S7_0.2.1