Title: | Clustering of MS2 Spectra for Metabolite Identification |
---|---|
Description: | CluMSID is a tool that aids the identification of features in untargeted LC-MS/MS analysis by the use of MS2 spectra similarity and unsupervised statistical methods. It offers functions for a complete and customisable workflow from raw data to visualisations and is interfaceable with the xmcs family of preprocessing packages. |
Authors: | Tobias Depke [aut, cre], Raimo Franke [ctb], Mark Broenstrup [ths] |
Maintainer: | Tobias Depke <[email protected]> |
License: | MIT + file LICENSE |
Version: | 1.23.0 |
Built: | 2024-11-30 03:23:20 UTC |
Source: | https://github.com/bioc/CluMSID |
MS2spectrum
and
pseudospectrum
objectsAccessor functions for individual slots of
MS2spectrum
and
pseudospectrum
objects
accessID(x) accessAnnotation(x) accessPrecursor(x) accessRT(x) accessPolarity(x) accessSpectrum(x) accessNeutralLosses(x)
accessID(x) accessAnnotation(x) accessPrecursor(x) accessRT(x) accessPolarity(x) accessSpectrum(x) accessNeutralLosses(x)
x |
An object of class |
The value of the respective slot of the object
(id
, annotation
, precursor
,
rt
, spectrum
, neutral_losses
)
load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) accessID(annotatedSpeclist[[1]]) load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) accessAnnotation(annotatedSpeclist[[1]]) load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) accessPrecursor(annotatedSpeclist[[1]]) load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) accessRT(annotatedSpeclist[[1]]) load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) accessPolarity(annotatedSpeclist[[1]]) load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) accessSpectrum(annotatedSpeclist[[1]]) load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) accessNeutralLosses(annotatedSpeclist[[1]])
load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) accessID(annotatedSpeclist[[1]]) load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) accessAnnotation(annotatedSpeclist[[1]]) load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) accessPrecursor(annotatedSpeclist[[1]]) load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) accessRT(annotatedSpeclist[[1]]) load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) accessPolarity(annotatedSpeclist[[1]]) load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) accessSpectrum(annotatedSpeclist[[1]]) load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) accessNeutralLosses(annotatedSpeclist[[1]])
MS2spectrum
objectsaddAnnotations
is used to add annotations that have been assigned
externally, e.g. by library search, to a list of MS2spectrum
objects
as produced by extractMS2spectra
and mergeSpecList
.
addAnnotations(featlist, annolist, annotationColumn = 4)
addAnnotations(featlist, annolist, annotationColumn = 4)
featlist |
A list of |
annolist |
A list of annotations, either as a |
annotationColumn |
The column of |
A list of MS2spectrum
objects as produced by
extractMS2spectra
and mergeSpecList
with external
annotations added to the annotation
slot of each MS2spectrum
object.
load(file = system.file("extdata", "featlist.RData", package = "CluMSIDdata")) addAnnotations(featlist, system.file("extdata", "post_anno.csv", package = "CluMSIDdata"), annotationColumn = 4)
load(file = system.file("extdata", "featlist.RData", package = "CluMSIDdata")) addAnnotations(featlist, system.file("extdata", "post_anno.csv", package = "CluMSIDdata"), annotationColumn = 4)
Convert spectra from MSnbase classes
as.MS2spectrum(x)
as.MS2spectrum(x)
x |
An object of class MS2spectrum
#Load a "Spectrum2" object from MSnbase library(MSnbase) sp <- itraqdata[["X1"]] #Convert this object to "MS2spectrum" class new_sp <- as.MS2spectrum(sp) #Or alternatively: new_sp <- as(sp, "MS2spectrum")
#Load a "Spectrum2" object from MSnbase library(MSnbase) sp <- itraqdata[["X1"]] #Convert this object to "MS2spectrum" class new_sp <- as.MS2spectrum(sp) #Or alternatively: new_sp <- as(sp, "MS2spectrum")
cossim()
calculates the cosine of the spectral constrast angle as a
measure for the similarity of two spectra.
cossim(x, y, type = c("spectrum", "neutral_losses"), mzTolerance = 1e-05) ## S4 method for signature 'MS2spectrum,MS2spectrum' cossim(x, y, type = c("spectrum", "neutral_losses"), mzTolerance = 1e-05) ## S4 method for signature 'pseudospectrum,pseudospectrum' cossim(x, y, type = c("spectrum", "neutral_losses"), mzTolerance = 1e-05)
cossim(x, y, type = c("spectrum", "neutral_losses"), mzTolerance = 1e-05) ## S4 method for signature 'MS2spectrum,MS2spectrum' cossim(x, y, type = c("spectrum", "neutral_losses"), mzTolerance = 1e-05) ## S4 method for signature 'pseudospectrum,pseudospectrum' cossim(x, y, type = c("spectrum", "neutral_losses"), mzTolerance = 1e-05)
x , y
|
MS2 spectra, either as |
type |
Whether similarity between spectra ( |
mzTolerance |
The m/z tolerance used for merging. If two fragment peaks
are within tolerance, they are regarded as the same. Defaults to
|
The cosine similarity of x
and y
x = MS2spectrum,y = MS2spectrum
: cossim
method for
MS2spectrum
objects
x = pseudospectrum,y = pseudospectrum
: cossim
method for
pseudospectrum
objects
load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) cossim(annotatedSpeclist[[1]], annotatedSpeclist[[2]])
load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) cossim(annotatedSpeclist[[1]], annotatedSpeclist[[2]])
distanceMatrix()
creates a distance matrix from a list of MS2
spectra, MS1 pseudospectra or neutral loss patterns by pairwise comparison
using the specified distance function. This distance matrix is the basis for
CluMSID's data mining functions.
distanceMatrix(speclist, distFun = "cossim", type = c("spectrum", "neutral_losses"), mz_tolerance = 1e-05)
distanceMatrix(speclist, distFun = "cossim", type = c("spectrum", "neutral_losses"), mz_tolerance = 1e-05)
speclist |
A list of |
distFun |
The distance function to be used. At the moment, only
|
type |
|
mz_tolerance |
The m/z tolerance to be used for merging, default
is |
A numeric length(speclist)
by length(speclist)
matrix
containing pairwise distances (1 - similarity) between all features in
speclist
. Row and column names are taken from the id
slot
or, if present, pasted from the id
and annotation
slots of
the MS2spectrum
or
pseudospectrum
objects.
load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) distanceMatrix(annotatedSpeclist[1:20])
load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) distanceMatrix(annotatedSpeclist[1:20])
extractMS2spectra()
is used to extract MS2 spectra from raw data
files, e.g. mzXML files.
extractMS2spectra(MSfile, min_peaks = 2, recalibrate_precursor = FALSE, RTlims = NULL)
extractMS2spectra(MSfile, min_peaks = 2, recalibrate_precursor = FALSE, RTlims = NULL)
MSfile |
An LC-MS/MS raw data file in one of the non-proprietary
formats that can be parsed by |
min_peaks |
Minimum number of peaks in MS2 spectrum, defaults to
|
recalibrate_precursor |
Logical, defaults to |
RTlims |
Retention time interval for the extraction of spectra. Provide
as numeric vector of length 2. Spectra with retention time <
|
A list
with objects of class MS2spectrum
, containing
MS2 spectra extracted from the raw data.
my_spectra <- extractMS2spectra(MSfile = system.file("extdata", "PoolA_R_SE.mzXML", package = "CluMSIDdata"), min_peaks = 4, RTlims = c(0,10))
my_spectra <- extractMS2spectra(MSfile = system.file("extdata", "PoolA_R_SE.mzXML", package = "CluMSIDdata"), min_peaks = 4, RTlims = c(0,10))
extractPseudospectra()
is used to extract MS1 pseudospectra from
CAMERA output.
extractPseudospectra(x, min_peaks = 1, intensity_columns = NULL)
extractPseudospectra(x, min_peaks = 1, intensity_columns = NULL)
x |
CAMERA output that contains information on pseudospectra.
Can either be of class |
min_peaks |
Minimum number of peaks in pseudospectrum, defaults to
|
intensity_columns |
Numeric, defaults to |
A list of pseudospectra, stored as objects of class
pseudospectrum
, analogous to the output of
extractMS2spectra
.
pstable <- readr::read_delim(file = system.file("extdata", "TD035_XCMS.annotated.diffreport.tsv", package = "CluMSIDdata"), delim = "\t") pseudospeclist <- extractPseudospectra(pstable, min_peaks = 2)
pstable <- readr::read_delim(file = system.file("extdata", "TD035_XCMS.annotated.diffreport.tsv", package = "CluMSIDdata"), delim = "\t") pseudospeclist <- extractPseudospectra(pstable, min_peaks = 2)
data.frame
with feature information from list of
MS2spectrum
objectsfeatureList
generates a data.frame
that contains feature ID,
precurosur m/z and retention time for all features contained in a
list of MS2spectrum
objects as produced by extractMS2spectra
and mergeSpecList
. featureList
is used internally by
writeFeaturelist
.
featureList(featlist)
featureList(featlist)
featlist |
A list of |
Although originally designed for lists of MS2spectrum
objects, the function also works with lists of pseudospectrum
objects. In this case, NA
is given for precursor m/z.
A data.frame
that contains feature ID, precurosur m/z
(if available) and retention time
load(file = system.file("extdata", "featlist.RData", package = "CluMSIDdata")) pre_anno <- featureList(featlist)
load(file = system.file("extdata", "featlist.RData", package = "CluMSIDdata")) pre_anno <- featureList(featlist)
findFragment
is used to find spectra that contain a specific fragment
ion. Its sister function is findNL
, which finds specific
neutral losses. Both functions work analogous to getSpectrum
.
findFragment(featlist, mz, tolerance = 1e-05)
findFragment(featlist, mz, tolerance = 1e-05)
featlist |
a list that contains only objects of class
|
mz |
The mass-to-charge ratio of the fragment ion of interest. |
tolerance |
The m/z tolerance for the fragment ion search.
Default is |
If the respective fragment is only found in one spectrum, the output
is an object of class MS2spectrum
; if it is found in
more than one spectrum, the output is a list of
MS2spectrum
objects.
load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) putativeAQs <- findFragment(annotatedSpeclist, 159.068)
load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) putativeAQs <- findFragment(annotatedSpeclist, 159.068)
findNL
is used to find spectra that contain a specific neutral loss.
Its sister function is findFragment
, which finds specific
fragment ions. Both functions work analogous to getSpectrum
.
findNL(featlist, mz, tolerance = 1e-05)
findNL(featlist, mz, tolerance = 1e-05)
featlist |
a list that contains only objects of class
|
mz |
The mass-to-charge ratio of the neutral loss of interest. |
tolerance |
The m/z tolerance for the neutral loss search.
Default is |
If the respective neutral loss is only found in one spectrum, the
output is an object of class MS2spectrum
; if it is
found in more than one spectrum, the output is a list of
MS2spectrum
objects.
load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) findNL(annotatedSpeclist, 212.009)
load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) findNL(annotatedSpeclist, 212.009)
getSimilarities
calculates the similarities of one spectrum or
neutral loss pattern to a set of other spectra or neutral loss patterns.
getSimilarities(spec, speclist, type = c("spectrum", "neutral_losses"), hits_only = FALSE)
getSimilarities(spec, speclist, type = c("spectrum", "neutral_losses"), hits_only = FALSE)
spec |
The spectrum to be compared to other spectra. Can be either an
object of class |
speclist |
The set of spectra to which |
type |
Specifies whether MS2 spectra or neutral loss patterns are to be
compared. Must be either |
hits_only |
Logical that indicates whether the result should contain only similarities greater than zero. |
A named vector with similarities of spec
to all spectra or
neutral loss patterns in speclist
.
load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) getSimilarities(annotatedSpeclist[[137]], annotatedSpeclist, hits_only = TRUE)
load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) getSimilarities(annotatedSpeclist[[137]], annotatedSpeclist, hits_only = TRUE)
As accessing S4 objects within lists is not trivial, getSpectrum
can
be used to access individual or several MS2spectrum
objects by their slot entries.
getSpectrum(featlist, slot, what, mz.tol = 1e-05, rt.tol = 30)
getSpectrum(featlist, slot, what, mz.tol = 1e-05, rt.tol = 30)
featlist |
a list that contains only objects of class
|
slot |
The slot to be searched (invalid slot arguments will produce errors). Possible values are:
|
what |
the search term or number, must be character for |
mz.tol |
the tolerance used for precursor ion *m/z* searches, defaults
to |
rt.tol |
the tolerance used for precursor ion retention time searches, defaults to 30s; high values can be used to specify retention time ranges (see vignette for example) |
If the only one spectrum matches the search criteria, the output is
an object of class MS2spectrum
; if more than one
spectrum matches, the output is a list of MS2spectrum
objects.
load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) getSpectrum(annotatedSpeclist, "annotation", "pyocyanin") getSpectrum(annotatedSpeclist, "id", "M244.17T796.4") getSpectrum(annotatedSpeclist, "precursor", 286.18, mz.tol = 1E-03) six_eight <- getSpectrum(annotatedSpeclist, "rt", 420, rt.tol = 60)
load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) getSpectrum(annotatedSpeclist, "annotation", "pyocyanin") getSpectrum(annotatedSpeclist, "id", "M244.17T796.4") getSpectrum(annotatedSpeclist, "precursor", 286.18, mz.tol = 1E-03) six_eight <- getSpectrum(annotatedSpeclist, "rt", 420, rt.tol = 60)
HCplot()
performs hierarchical clustering
of spectral similarity data using average linkage
as agglomeration criterion like HCtbl
and generates either a circular dendrogram or a
combination of dendrogram and heatmap.
HCplot(distmat, h = 0.95, type = c("dendrogram", "heatmap"), ...)
HCplot(distmat, h = 0.95, type = c("dendrogram", "heatmap"), ...)
distmat |
A distance matrix as generated by
|
h |
Height where the tree is to be cut, defaults to |
type |
Specifies which visualisation is to be generated:
|
... |
Additional graphical parameters passed to
|
A plot as specified by type
.
load(file = system.file("extdata", "distmat.RData", package = "CluMSIDdata")) HCplot(distmat[1:50,1:50], h = 0.8, type = "heatmap")
load(file = system.file("extdata", "distmat.RData", package = "CluMSIDdata")) HCplot(distmat[1:50,1:50], h = 0.8, type = "heatmap")
HCtbl()
performs hierarchical clustering
of spectral similarity data using average linkage
as agglomeration criterion.
HCtbl(distmat, h = 0.95)
HCtbl(distmat, h = 0.95)
distmat |
A distance matrix as generated by
|
h |
Height where the tree is to be cut, defaults to |
A data.frame
with name and cluster ID for each
feature in distmat
.
load(file = system.file("extdata", "distmat.RData", package = "CluMSIDdata")) my_HCtbl <- HCtbl(distmat[1:50,1:50], h = 0.8)
load(file = system.file("extdata", "distmat.RData", package = "CluMSIDdata")) my_HCtbl <- HCtbl(distmat[1:50,1:50], h = 0.8)
MDSplot()
is used to generate multidimensional scaling plots from
spectral similarity data. An interactive visualisation can be produced using
plotly.
MDSplot(distmat, interactive = FALSE, highlight_annotated = FALSE, ...)
MDSplot(distmat, interactive = FALSE, highlight_annotated = FALSE, ...)
distmat |
A distance matrix as generated by
|
interactive |
Logical, defaults to |
highlight_annotated |
Logical, defaults to |
... |
Additional arguments passed to |
An MDS plot generated with the help of
cmdscale
,
ggplot
and, if interactive,
ggplotly
.
load(file = system.file("extdata", "distmat.RData", package = "CluMSIDdata")) MDSplot(distmat, highlight_annotated = TRUE)
load(file = system.file("extdata", "distmat.RData", package = "CluMSIDdata")) MDSplot(distmat, highlight_annotated = TRUE)
mergeMS2spectra
is used to merge MS2 spectra that come from the same
precursor. It does so either by grouping spectra of the same precursor
m/z that fall into a defined retention time window
(rt_tolerance
) or by grouping spectra to peaks from an externally
supplied peak table. Please see the vignette for more details.
mergeMS2spectra(ms2list, mz_tolerance = 1e-05, rt_tolerance = 30, peaktable = NULL, exclude_unmatched = FALSE)
mergeMS2spectra(ms2list, mz_tolerance = 1e-05, rt_tolerance = 30, peaktable = NULL, exclude_unmatched = FALSE)
ms2list |
A |
mz_tolerance |
The m/z tolerance to be used for merging, default
is |
rt_tolerance |
The retention time tolerance used for merging features.
If used without a peak table, |
peaktable |
An external peak table, e.g. from XCMS, that serves as a
template for grouping spectra. The peaktable must be a three-column
|
exclude_unmatched |
If an external peak table is used: Should spectra that do not match to any peak/feature in the peak table be exclude from the resulting list? |
A merged list of MS2spectrum
objects.
my_spectra <- extractMS2spectra(MSfile = system.file("extdata", "PoolA_R_SE.mzXML", package = "CluMSIDdata"), min_peaks = 4, RTlims = c(0,5)) my_merged_spectra <- mergeMS2spectra(my_spectra, rt_tolerance = 20)
my_spectra <- extractMS2spectra(MSfile = system.file("extdata", "PoolA_R_SE.mzXML", package = "CluMSIDdata"), min_peaks = 4, RTlims = c(0,5)) my_merged_spectra <- mergeMS2spectra(my_spectra, rt_tolerance = 20)
A custom S4 class for MS2 spectra, neutral loss patterns and respective metainformation
## S4 method for signature 'MS2spectrum' show(object) ## S4 method for signature 'MS2spectrum' precursorMz(object) ## S4 method for signature 'MS2spectrum' rtime(object) ## S4 method for signature 'MS2spectrum' intensity(object) ## S4 method for signature 'MS2spectrum' mz(object) ## S4 method for signature 'MS2spectrum,ANY' peaksCount(object)
## S4 method for signature 'MS2spectrum' show(object) ## S4 method for signature 'MS2spectrum' precursorMz(object) ## S4 method for signature 'MS2spectrum' rtime(object) ## S4 method for signature 'MS2spectrum' intensity(object) ## S4 method for signature 'MS2spectrum' mz(object) ## S4 method for signature 'MS2spectrum,ANY' peaksCount(object)
object |
An object of class |
Prints information from the object slots with exception of 'spectrum' and 'neutral_losses' where only a summary is given.
show
: A show generic for MS2spectra
.
precursorMz
: Method forMSnbase::precursorMz
for MS2spectrum
objects. Accesses precursor
slot and returns precursor m/z as a numeric.
rtime
: Method forMSnbase::rtime
for MS2spectrum
objects. Accesses rt
slot and returns retention time as a numeric.
intensity
: Method forMSnbase::intensity
for MS2spectrum
objects. Accesses spectrum
slot and returns the intensity column as a numeric vector.
mz
: Method forMSnbase::mz
for MS2spectrum
objects. Accesses spectrum
slot and returns the m/z column as a numeric vector.
peaksCount
: Method forMSnbase::mz
for MS2spectrum
objects. Accesses spectrum
slot and returns the number of peaks as a numeric.
id
a character string similar to the ID used by XCMSonline or the ID given in a predefined peak list
annotation
a character string containing a user-defined annotation, defaults to empty
precursor
(median) m/z of the spectrum's precursor ion
rt
(median) retention time of the spectrum's precursor ion
polarity
the ionisation polarity,
"positive"
or "negative"
spectrum
the actual MS2 spectrum as two-column matrix (column 1 is (median) m/z, column 2 is (median) intensity of the product ions)
neutral_losses
a neutral loss pattern generated by subtracting the
product ion mass-to-charge ratios from the precursor m/z in a
matrix format analogous to the spectrum
slot
networkplot()
is used to generate correlation networks from
spectral similarity data. An interactive visualisation can be produced using
plotly.
networkplot(distmat, interactive = FALSE, show_labels = FALSE, label_size = 1.5, highlight_annotated = FALSE, min_similarity = 0.1, exclude_singletons = FALSE)
networkplot(distmat, interactive = FALSE, show_labels = FALSE, label_size = 1.5, highlight_annotated = FALSE, min_similarity = 0.1, exclude_singletons = FALSE)
distmat |
A distance matrix as generated by
|
interactive |
Logical, defaults to |
show_labels |
Logical, defaults to |
label_size |
Numeric, defaults to |
highlight_annotated |
Logical, defaults to |
min_similarity |
Numeric, defaults to |
exclude_singletons |
Logical, defaults to |
A network plot generated with the help of
network
, ggnet2
and, if
interactive, ggplotly
. Edge weights correspond to
spectral similarities.
load(file = system.file("extdata", "distmat.RData", package = "CluMSIDdata")) networkplot(distmat[1:50,1:50], show_labels = TRUE, exclude_singletons = TRUE)
load(file = system.file("extdata", "distmat.RData", package = "CluMSIDdata")) networkplot(distmat[1:50,1:50], show_labels = TRUE, exclude_singletons = TRUE)
OPTICSplot()
performs density-based clustering
of spectral similarity data using the OPTICS algorithm like
OPTICStbl
and creates a reachability distance plot.
OPTICSplot(distmat, eps = 10000, minPts = 3, eps_cl = 0.5, ...)
OPTICSplot(distmat, eps = 10000, minPts = 3, eps_cl = 0.5, ...)
distmat |
A distance matrix as generated by
|
eps |
OPTICS parameters, see |
minPts |
OPTICS parameters, see |
eps_cl |
The reachability distance used for cluster determination,
see |
... |
Additional graphical parameters to be passed to |
The function internally uses optics
and extractDBSCAN
from the dbscan
package.
A reachability distance plot as visualisation of OPTICS clustering, see codeextractDBSCAN.
OPTICStbl
load(file = system.file("extdata", "distmat.RData", package = "CluMSIDdata")) OPTICSplot(distmat[1:50,1:50], eps_cl = 0.7)
load(file = system.file("extdata", "distmat.RData", package = "CluMSIDdata")) OPTICSplot(distmat[1:50,1:50], eps_cl = 0.7)
OPTICStbl()
performs density-based clustering
of spectral similarity data using the OPTICS algorithm.
OPTICStbl(distmat, eps = 10000, minPts = 3, eps_cl = 0.5)
OPTICStbl(distmat, eps = 10000, minPts = 3, eps_cl = 0.5)
distmat |
A distance matrix as generated by
|
eps , minPts
|
OPTICS parameters, see |
eps_cl |
The reachability distance used for cluster determination,
see |
The function internally uses optics
and extractDBSCAN
from the dbscan
package.
A data.frame
with feature name, cluster ID and
OPTICS order for each feature in distmat
.
OPTICSplot
load(file = system.file("extdata", "distmat.RData", package = "CluMSIDdata")) my_OTPICStbl <- OPTICStbl(distmat[1:50,1:50], eps_cl = 0.7)
load(file = system.file("extdata", "distmat.RData", package = "CluMSIDdata")) my_OTPICStbl <- OPTICStbl(distmat[1:50,1:50], eps_cl = 0.7)
A custom S4 class for MS1 pseudospectra and respective metainformation
id
a the "pcgroup"
number assigned by CAMERA
annotation
a character string containing a user-defined annotation, defaults to empty
rt
(median) retention time of the ions contained in the pseudospectrum
spectrum
the actual MS1 pseudospectrum as two-column matrix (column 1 is (median) m/z, column 2 is (median) intensity of the ions)
specplot
creates a very basic plot of MS2 spectra from
MS2spectrum
or pseudospectrum
objects.
specplot(spec, ...)
specplot(spec, ...)
spec |
An object of class |
... |
Additional graphical parameters to be passed to |
A plot of the MS2 spectrum saved in the spectrum
slot of
spec
.
load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) specplot(annotatedSpeclist[[1]])
load(file = system.file("extdata", "annotatedSpeclist.RData", package = "CluMSIDdata")) specplot(annotatedSpeclist[[1]])
Using splitPolarities
, spectra with different polarities from the same
run can be separated, e.g. when processing spectra recorded with
polarity-switching.
splitPolarities(ms2list, polarity = c("positive", "negative"))
splitPolarities(ms2list, polarity = c("positive", "negative"))
ms2list |
A list of |
polarity |
The polarity of spectra to be analysed, must be
|
A list of MS2spectrum
objects that contains only
spectra with the given polarity
.
my_spectra <- extractMS2spectra(MSfile = system.file("extdata", "PoolA_R_SE.mzXML", package = "CluMSIDdata"), min_peaks = 4, RTlims = c(0,5)) my_positive_spectra <- splitPolarities(my_spectra, "positive")
my_spectra <- extractMS2spectra(MSfile = system.file("extdata", "PoolA_R_SE.mzXML", package = "CluMSIDdata"), min_peaks = 4, RTlims = c(0,5)) my_positive_spectra <- splitPolarities(my_spectra, "positive")
MS2spectrum
objectswriteFeaturelist
uses featureList
to generate a
data.frame
that contains feature ID, precurosur m/z and
retention time for all features contained in a list of MS2spectrum
objects as produced by extractMS2spectra
and mergeSpecList
and
writes it to a csv file.
writeFeaturelist(featlist, filename = "pre_anno.csv")
writeFeaturelist(featlist, filename = "pre_anno.csv")
featlist |
A list of |
filename |
The desired file name of the csv file, default is
|
Although originally designed for lists of MS2spectrum
objects, the function also works with lists of pseudospectrum
objects. In this case, NA
is given for precursor m/z.
A csv file that contains feature ID, precurosur m/z and
retention time. The file has a header but no row names and is separated by
','
.
load(file = system.file("extdata", "featlist.RData", package = "CluMSIDdata")) writeFeaturelist(featlist, filename = "pre_anno.csv")
load(file = system.file("extdata", "featlist.RData", package = "CluMSIDdata")) writeFeaturelist(featlist, filename = "pre_anno.csv")