peakPantheR 1.21.0
Package: peakPantheR
Authors: Arnaud Wolfer
Package for Peak Picking and ANnoTation of High resolution Experiments in R,
implemented in R
and Shiny
peakPantheR
implements functions to detect, integrate and report pre-defined
features in MS files (e.g. compounds, fragments, adducts, …).
It is designed for:
multiple
compounds in one
file at a timemultiple
compounds in multiple
files in parallel
, store
results in a single
objectpeakPantheR
can process LC/MS data files in NetCDF, mzML/mzXML and
mzData format as data import is achieved using Bioconductor’s
mzR package.
To install peakPantheR
from Bioconductor:
if (!requireNamespace("BiocManager", quietly = TRUE))
install.packages("BiocManager")
BiocManager::install("peakPantheR")
Install the development version of peakPantheR
directly from GitHub with:
# Install devtools
if(!require("devtools")) install.packages("devtools")
devtools::install_github("phenomecentre/peakPantheR")
To get started peakPantheR
’s graphical user interface
implements all the functions to detect and integrate multiple compounds in
multiple files in parallel and store results in a single object. It
can be employed to integrate a large number of expected features across a
dataset:
library(peakPantheR)
peakPantheR_start_GUI(browser = TRUE)
# To exit press ESC in the command line
The GUI is to be preferred to understand the methodology, select the best parameters on a subset of the samples before running the command line, or to visually explore results.
If a very high number of samples and features is to be processed,
peakpantheR
’s command line functions are more efficient, as they can
be integrated in scripts and the reporting automated.
Both real time and parallel compound integration require a common set of information:
netCDF
/ mzML
MS file(s)RT
/ m/z
window) for each compound.For demonstration purpose we can annotate a set a set of raw MS spectra (in NetCDF format) provided by the faahKO package. Briefly, this subset of the data from (Saghatelian et al. 2004) invesigate the metabolic consequences of knocking out the fatty acid amide hydrolase (FAAH) gene in mice. The dataset consists of samples from the spinal cords of 6 knock-out and 6 wild-type mice. Each file contains data in centroid mode acquired in positive ion mode form 200-600 m/z and 2500-4500 seconds.
Below we install the faahKO package and locate raw CDF files of interest:
if (!requireNamespace("BiocManager", quietly = TRUE))
install.packages("BiocManager")
BiocManager::install("faahKO")
library(faahKO)
## file paths
input_spectraPaths <- c(system.file('cdf/KO/ko15.CDF', package = "faahKO"),
system.file('cdf/KO/ko16.CDF', package = "faahKO"),
system.file('cdf/KO/ko18.CDF', package = "faahKO"))
input_spectraPaths
#> [1] "/home/biocbuild/bbs-3.21-bioc/R/site-library/faahKO/cdf/KO/ko15.CDF"
#> [2] "/home/biocbuild/bbs-3.21-bioc/R/site-library/faahKO/cdf/KO/ko16.CDF"
#> [3] "/home/biocbuild/bbs-3.21-bioc/R/site-library/faahKO/cdf/KO/ko18.CDF"
Expected regions of interest (targeted features) are specified using the following information:
cpdID
(character)cpdName
(character)rtMin
(sec)rtMax
(sec)rt
(sec, optional / NA
)mzMin
(m/z)mzMax
(m/z)mz
(m/z, optional / NA
)Below we define 2 features of interest that are present in the faahKO dataset and can be employed in subsequent vignettes:
# targetFeatTable
input_targetFeatTable <- data.frame(matrix(vector(), 2, 8, dimnames=list(c(),
c("cpdID", "cpdName", "rtMin", "rt", "rtMax", "mzMin",
"mz", "mzMax"))), stringsAsFactors=FALSE)
input_targetFeatTable[1,] <- c("ID-1", "Cpd 1", 3310., 3344.888, 3390., 522.194778,
522.2, 522.205222)
input_targetFeatTable[2,] <- c("ID-2", "Cpd 2", 3280., 3385.577, 3440., 496.195038,
496.2, 496.204962)
input_targetFeatTable[,c(1,3:8)] <- sapply(input_targetFeatTable[,c(1,3:8)],
as.numeric)
#> Warning in lapply(X = X, FUN = FUN, ...): NAs introduced by coercion
cpdID | cpdName | rtMin | rt | rtMax | mzMin | mz | mzMax |
---|---|---|---|---|---|---|---|
NA | Cpd 1 | 3310 | 3344.888 | 3390 | 522.194778 | 522.2 | 522.205222 |
NA | Cpd 2 | 3280 | 3385.577 | 3440 | 496.195038 | 496.2 | 496.204962 |
While the command line functions accept Data.Frame and vectors as input, the
graphical user interface (GUI) will read the same information from a set of
.csv
files, or an already set-up peakPantheRAnnotation
object in .RData
format.
We can now generate GUI input files for the faahKO example dataset presented previously:
A peakPantheRAnnotation
(previously annotated or not) can be passed as input
in a .RData
file. The peakPantheRAnnotation
object must be named
annotationObject:
library(faahKO)
# Define the file paths (3 samples)
input_spectraPaths <- c(system.file('cdf/KO/ko15.CDF', package = "faahKO"),
system.file('cdf/KO/ko16.CDF', package = "faahKO"),
system.file('cdf/KO/ko18.CDF', package = "faahKO"))
# Define the targeted features (2 features)
input_targetFeatTable <- data.frame(matrix(vector(), 2, 8, dimnames=list(c(),
c("cpdID", "cpdName", "rtMin", "rt", "rtMax", "mzMin",
"mz", "mzMax"))), stringsAsFactors=FALSE)
input_targetFeatTable[1,] <- c("ID-1", "Cpd 1", 3310., 3344.888, 3390.,
522.194778, 522.2, 522.205222)
input_targetFeatTable[2,] <- c("ID-2", "Cpd 2", 3280., 3385.577, 3440.,
496.195038, 496.2, 496.204962)
input_targetFeatTable[,3:8] <- sapply(input_targetFeatTable[,3:8], as.numeric)
# Define some random compound and spectra metadata
# cpdMetadata
input_cpdMetadata <- data.frame(matrix(data=c('a','b',1,2), nrow=2, ncol=2,
dimnames=list(c(), c('testcol1','testcol2')),
byrow=FALSE), stringsAsFactors=FALSE)
# spectraMetadata
input_spectraMetadata <- data.frame(matrix(data=c('c','d','e',3,4,5), nrow=3,
ncol=2,
dimnames=list(c(),c('testcol1','testcol2')),
byrow=FALSE), stringsAsFactors=FALSE)
# Initialise a simple peakPantheRAnnotation object
# [3 files, 2 features, no uROI, no FIR]
initAnnotation <- peakPantheRAnnotation(spectraPaths=input_spectraPaths,
targetFeatTable=input_targetFeatTable,
cpdMetadata=input_cpdMetadata,
spectraMetadata=input_spectraMetadata)
# Rename and save the annotation to disk
annotationObject <- initAnnotation
save(annotationObject,
file = './example_annotation_ppR_UI.RData',
compress=TRUE)
Another input option for the GUI input consists of a set of .csv
files.
Targeted features are defined in a .csv
with as rows each feature to target
(the first row must be the column name), and as columns the fit parameters to
use. At minimum the following parameters must be defined:
cpdID
, cpdName
, rtMin
, rt
, rtMax
, mzMin
, mz
, mzMax
If uROI
and FIR
are to be set, the following columns must be provided:
cpdID
, cpdName
, ROI_rt
, ROI_mz
, ROI_rtMin
, ROI_rtMax
, ROI_mzMin
,
ROI_mzMax
, uROI_rtMin
, uROI_rtMax
, uROI_mzMin
, uROI_mzMax
, uROI_rt
,
uROI_mz
, FIR_rtMin
, FIR_rtMax
, FIR_mzMin
, FIR_mzMax
# Define targeted features without uROI and FIR (2 features)
input_targetFeatTable <- data.frame(matrix(vector(), 2, 8, dimnames=list(c(),
c("cpdID", "cpdName", "rtMin", "rt", "rtMax", "mzMin",
"mz", "mzMax"))), stringsAsFactors=FALSE)
input_targetFeatTable[1,] <- c("ID-1", "Cpd 1", 3310., 3344.888, 3390.,
522.194778, 522.2, 522.205222)
input_targetFeatTable[2,] <- c("ID-2", "Cpd 2", 3280., 3385.577, 3440.,
496.195038, 496.2, 496.204962)
input_targetFeatTable[,3:8] <- sapply(input_targetFeatTable[,3:8], as.numeric)
# save to disk
write.csv(input_targetFeatTable,
file = './1-fitParams_example_UI.csv',
row.names = FALSE)
cpdID | cpdName | rtMin | rt | rtMax | mzMin | mz | mzMax |
---|---|---|---|---|---|---|---|
ID-1 | Cpd 1 | 3310 | 3344.888 | 3390 | 522.194778 | 522.2 | 522.205222 |
ID-2 | Cpd 2 | 3280 | 3385.577 | 3440 | 496.195038 | 496.2 | 496.204962 |
It is possible to select the files on disk directly through the GUI, or to
select a .csv
file containing each file path as well as spectra metadata.
Each row correspond to a different spectra (the first row must define the
column names) while columns correspond to the path on disk and metadata. At
minimum a column filepath
must be present, with subsequent columns defining
metadata properties.
# Define the spectra paths and metada
input_spectraMeta <- data.frame(matrix(vector(), 3, 3,
dimnames=list(c(),c("filepath","testcol1","testcol2"))),
stringsAsFactors=FALSE)
input_spectraMeta[1,] <- c(system.file('cdf/KO/ko15.CDF', package = "faahKO"),
"c", 3)
input_spectraMeta[2,] <- c(system.file('cdf/KO/ko16.CDF', package = "faahKO"),
"d", 4)
input_spectraMeta[3,] <- c(system.file('cdf/KO/ko18.CDF', package = "faahKO"),
"e", 5)
# save to disk
write.csv(input_spectraMeta,
file = './2-spectraMetaWPath_example_UI.csv',
row.names = FALSE)
filepath | testcol1 |
---|---|
/home/biocbuild/bbs-3.21-bioc/R/site-library/faahKO/cdf/KO/ko15.CDF | c |
/home/biocbuild/bbs-3.21-bioc/R/site-library/faahKO/cdf/KO/ko16.CDF | d |
/home/biocbuild/bbs-3.21-bioc/R/site-library/faahKO/cdf/KO/ko18.CDF | e |
testcol2 |
---|
3 |
4 |
5 |
It is possible to define some feature metadata, with targeted features as rows (same order as the fitting parameters, first row defining the column names), and as columns the metadata.
# Define the feature metada
input_featMeta <- data.frame(matrix(vector(), 2, 2,
dimnames=list(c(),c("testcol1","testcol2"))),
stringsAsFactors=FALSE)
input_featMeta[1,] <- c("a", 1)
input_featMeta[2,] <- c("b", 2)
# save to disk
write.csv(input_featMeta,
file = './3-featMeta_example_UI.csv',
row.names = FALSE)
testcol1 | testcol2 |
---|---|
a | 1 |
b | 2 |
#> ─ Session info ───────────────────────────────────────────────────────────────
#> setting value
#> version R Under development (unstable) (2024-10-21 r87258)
#> os Ubuntu 24.04.1 LTS
#> system x86_64, linux-gnu
#> ui X11
#> language (EN)
#> collate C
#> ctype en_US.UTF-8
#> tz America/New_York
#> date 2024-10-29
#> pandoc 3.1.3 @ /usr/bin/ (via rmarkdown)
#>
#> ─ Packages ───────────────────────────────────────────────────────────────────
#> package * version date (UTC) lib source
#> abind 1.4-8 2024-09-12 [2] CRAN (R 4.5.0)
#> affy 1.85.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> affyio 1.77.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> AnnotationFilter 1.31.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> Biobase 2.67.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> BiocGenerics 0.53.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> BiocManager 1.30.25 2024-08-28 [2] CRAN (R 4.5.0)
#> BiocParallel * 1.41.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> BiocStyle * 2.35.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> bookdown 0.41 2024-10-16 [2] CRAN (R 4.5.0)
#> bslib 0.8.0 2024-07-29 [2] CRAN (R 4.5.0)
#> cachem 1.1.0 2024-05-16 [2] CRAN (R 4.5.0)
#> cli 3.6.3 2024-06-21 [2] CRAN (R 4.5.0)
#> clue 0.3-65 2023-09-23 [2] CRAN (R 4.5.0)
#> cluster 2.1.6 2023-12-01 [3] CRAN (R 4.5.0)
#> codetools 0.2-20 2024-03-31 [3] CRAN (R 4.5.0)
#> colorspace 2.1-1 2024-07-26 [2] CRAN (R 4.5.0)
#> crayon 1.5.3 2024-06-20 [2] CRAN (R 4.5.0)
#> DBI 1.2.3 2024-06-02 [2] CRAN (R 4.5.0)
#> DelayedArray 0.33.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> devtools 2.4.5 2022-10-11 [2] CRAN (R 4.5.0)
#> digest 0.6.37 2024-08-19 [2] CRAN (R 4.5.0)
#> doParallel 1.0.17 2022-02-07 [2] CRAN (R 4.5.0)
#> dplyr 1.1.4 2023-11-17 [2] CRAN (R 4.5.0)
#> DT 0.33 2024-04-04 [2] CRAN (R 4.5.0)
#> ellipsis 0.3.2 2021-04-29 [2] CRAN (R 4.5.0)
#> evaluate 1.0.1 2024-10-10 [2] CRAN (R 4.5.0)
#> faahKO * 1.45.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> fansi 1.0.6 2023-12-08 [2] CRAN (R 4.5.0)
#> fastmap 1.2.0 2024-05-15 [2] CRAN (R 4.5.0)
#> foreach 1.5.2 2022-02-02 [2] CRAN (R 4.5.0)
#> fs 1.6.4 2024-04-25 [2] CRAN (R 4.5.0)
#> generics 0.1.3 2022-07-05 [2] CRAN (R 4.5.0)
#> GenomeInfoDb 1.43.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> GenomeInfoDbData 1.2.13 2024-10-23 [2] Bioconductor
#> GenomicRanges 1.59.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> ggplot2 3.5.1 2024-04-23 [2] CRAN (R 4.5.0)
#> glue 1.8.0 2024-09-30 [2] CRAN (R 4.5.0)
#> gtable 0.3.6 2024-10-25 [2] CRAN (R 4.5.0)
#> highr 0.11 2024-05-26 [2] CRAN (R 4.5.0)
#> hms 1.1.3 2023-03-21 [2] CRAN (R 4.5.0)
#> htmltools 0.5.8.1 2024-04-04 [2] CRAN (R 4.5.0)
#> htmlwidgets 1.6.4 2023-12-06 [2] CRAN (R 4.5.0)
#> httpuv 1.6.15 2024-03-26 [2] CRAN (R 4.5.0)
#> httr 1.4.7 2023-08-15 [2] CRAN (R 4.5.0)
#> igraph 2.1.1 2024-10-19 [2] CRAN (R 4.5.0)
#> impute 1.81.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> IRanges 2.41.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> iterators 1.0.14 2022-02-05 [2] CRAN (R 4.5.0)
#> jquerylib 0.1.4 2021-04-26 [2] CRAN (R 4.5.0)
#> jsonlite 1.8.9 2024-09-20 [2] CRAN (R 4.5.0)
#> knitr 1.48 2024-07-07 [2] CRAN (R 4.5.0)
#> later 1.3.2 2023-12-06 [2] CRAN (R 4.5.0)
#> lattice 0.22-6 2024-03-20 [3] CRAN (R 4.5.0)
#> lazyeval 0.2.2 2019-03-15 [2] CRAN (R 4.5.0)
#> lifecycle 1.0.4 2023-11-07 [2] CRAN (R 4.5.0)
#> limma 3.63.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> magrittr 2.0.3 2022-03-30 [2] CRAN (R 4.5.0)
#> MALDIquant 1.22.3 2024-08-19 [2] CRAN (R 4.5.0)
#> MASS 7.3-61 2024-06-13 [3] CRAN (R 4.5.0)
#> MassSpecWavelet 1.73.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> Matrix 1.7-1 2024-10-18 [3] CRAN (R 4.5.0)
#> MatrixGenerics 1.19.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> matrixStats 1.4.1 2024-09-08 [2] CRAN (R 4.5.0)
#> memoise 2.0.1 2021-11-26 [2] CRAN (R 4.5.0)
#> MetaboCoreUtils 1.15.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> mime 0.12 2021-09-28 [2] CRAN (R 4.5.0)
#> miniUI 0.1.1.1 2018-05-18 [2] CRAN (R 4.5.0)
#> MsCoreUtils 1.19.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> MsExperiment 1.9.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> MsFeatures 1.15.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> MSnbase 2.33.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> MultiAssayExperiment 1.33.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> munsell 0.5.1 2024-04-01 [2] CRAN (R 4.5.0)
#> mzID 1.45.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> mzR 2.41.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> ncdf4 1.23 2024-08-17 [2] CRAN (R 4.5.0)
#> pander * 0.6.5 2022-03-18 [2] CRAN (R 4.5.0)
#> pcaMethods 1.99.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> peakPantheR * 1.21.0 2024-10-29 [1] Bioconductor 3.21 (R 4.5.0)
#> pillar 1.9.0 2023-03-22 [2] CRAN (R 4.5.0)
#> pkgbuild 1.4.5 2024-10-28 [2] CRAN (R 4.5.0)
#> pkgconfig 2.0.3 2019-09-22 [2] CRAN (R 4.5.0)
#> pkgload 1.4.0 2024-06-28 [2] CRAN (R 4.5.0)
#> plyr 1.8.9 2023-10-02 [2] CRAN (R 4.5.0)
#> preprocessCore 1.69.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> prettyunits 1.2.0 2023-09-24 [2] CRAN (R 4.5.0)
#> profvis 0.4.0 2024-09-20 [2] CRAN (R 4.5.0)
#> progress 1.2.3 2023-12-06 [2] CRAN (R 4.5.0)
#> promises 1.3.0 2024-04-05 [2] CRAN (R 4.5.0)
#> ProtGenerics 1.39.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> PSMatch 1.11.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> purrr 1.0.2 2023-08-10 [2] CRAN (R 4.5.0)
#> QFeatures 1.17.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> R6 2.5.1 2021-08-19 [2] CRAN (R 4.5.0)
#> RColorBrewer 1.1-3 2022-04-03 [2] CRAN (R 4.5.0)
#> Rcpp 1.0.13 2024-07-17 [2] CRAN (R 4.5.0)
#> remotes 2.5.0 2024-03-17 [2] CRAN (R 4.5.0)
#> reshape2 1.4.4 2020-04-09 [2] CRAN (R 4.5.0)
#> rlang 1.1.4 2024-06-04 [2] CRAN (R 4.5.0)
#> rmarkdown 2.28 2024-08-17 [2] CRAN (R 4.5.0)
#> S4Arrays 1.7.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> S4Vectors 0.45.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> sass 0.4.9 2024-03-15 [2] CRAN (R 4.5.0)
#> scales 1.3.0 2023-11-28 [2] CRAN (R 4.5.0)
#> sessioninfo 1.2.2 2021-12-06 [2] CRAN (R 4.5.0)
#> shiny 1.9.1 2024-08-01 [2] CRAN (R 4.5.0)
#> shinycssloaders 1.1.0 2024-07-30 [2] CRAN (R 4.5.0)
#> SparseArray 1.7.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> Spectra 1.17.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> statmod 1.5.0 2023-01-06 [2] CRAN (R 4.5.0)
#> stringi 1.8.4 2024-05-06 [2] CRAN (R 4.5.0)
#> stringr 1.5.1 2023-11-14 [2] CRAN (R 4.5.0)
#> SummarizedExperiment 1.37.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> svglite 2.1.3 2023-12-08 [2] CRAN (R 4.5.0)
#> systemfonts 1.1.0 2024-05-15 [2] CRAN (R 4.5.0)
#> tibble 3.2.1 2023-03-20 [2] CRAN (R 4.5.0)
#> tidyr 1.3.1 2024-01-24 [2] CRAN (R 4.5.0)
#> tidyselect 1.2.1 2024-03-11 [2] CRAN (R 4.5.0)
#> UCSC.utils 1.3.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> urlchecker 1.0.1 2021-11-30 [2] CRAN (R 4.5.0)
#> usethis 3.0.0 2024-07-29 [2] CRAN (R 4.5.0)
#> utf8 1.2.4 2023-10-22 [2] CRAN (R 4.5.0)
#> vctrs 0.6.5 2023-12-01 [2] CRAN (R 4.5.0)
#> vsn 3.75.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> xcms * 4.5.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> xfun 0.48 2024-10-03 [2] CRAN (R 4.5.0)
#> XML 3.99-0.17 2024-06-25 [2] CRAN (R 4.5.0)
#> xtable 1.8-4 2019-04-21 [2] CRAN (R 4.5.0)
#> XVector 0.47.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#> yaml 2.3.10 2024-07-26 [2] CRAN (R 4.5.0)
#> zlibbioc 1.53.0 2024-10-29 [2] Bioconductor 3.21 (R 4.5.0)
#>
#> [1] /tmp/RtmpUp8UHo/Rinst16108273a535d7
#> [2] /home/biocbuild/bbs-3.21-bioc/R/site-library
#> [3] /home/biocbuild/bbs-3.21-bioc/R/library
#>
#> ──────────────────────────────────────────────────────────────────────────────