Title: | Navigating mass spectral similarity in high-resolution MS/MS metabolomics data metabolomics data |
---|---|
Description: | MetCirc comprises a workflow to interactively explore high-resolution MS/MS metabolomics data. MetCirc uses the Spectra object infrastructure defined in the package Spectra that stores MS/MS spectra. MetCirc offers functionality to calculate similarity between precursors based on the normalised dot product, neutral losses or user-defined functions and visualise similarities in a circular layout. Within the interactive framework the user can annotate MS/MS features based on their similarity to (known) related MS/MS features. |
Authors: | Thomas Naake <[email protected]>, Johannes Rainer <[email protected]> and Emmanuel Gaquerel <[email protected]> |
Maintainer: | Thomas Naake <[email protected]> |
License: | GPL (>= 3) |
Version: | 1.37.0 |
Built: | 2024-12-09 06:16:27 UTC |
Source: | https://github.com/bioc/MetCirc |
cart2Polar
calculates polar coordinates from cartesian coordinates.
cart2Polar(x, y)
cart2Polar(x, y)
x |
|
y |
|
cart2Polar
is employed to translate cartesian coordinates
into polar coordinates especially in interactive shiny applications when
using hovering and clicking features.
cart2Polar
returns a list
of colar coordinates r and theta
Thomas Naake, [email protected]
x <- 1; y <- 1 cart2Polar(x, y)
x <- 1; y <- 1 cart2Polar(x, y)
circosLegend
plots a legend for circos plot using group names.
circosLegend(groupname, highlight = TRUE, colour = NULL, cex = 1)
circosLegend(groupname, highlight = TRUE, colour = NULL, cex = 1)
groupname |
|
highlight |
|
colour |
|
cex |
|
Internal use in shinyCircos
or outside of shinyCircos
to reproduce figures.
The function will open a new plot and display colours together with labels.
Thomas Naake, [email protected]
data("spectra", package = "MetCirc") ## create similarity matrix similarityMat <- Spectra::compareSpectra(sps_tissue[1:10], FUN = MsCoreUtils::ndotproduct, ppm = 20, m = 0.5, n = 2) rownames(similarityMat) <- colnames(similarityMat) <- sps_tissue$name[1:10] linkDf <- createLinkDf(similarityMatrix = similarityMat, sps = sps_tissue[1:10], condition = c("SPL", "LIM", "ANT", "STY"), lower = 0.01, upper = 1) ## cut link data.frame (here: only display links between groups) linkDf_cut <- cutLinkDf(linkDf, type = "inter") groupname <- c(as.character(linkDf_cut[, "spectrum1"]), as.character(linkDf_cut[, "spectrum2"])) groupname <- unique(groupname) ## plot legend circosLegend(groupname, highlight = TRUE, colour = NULL, cex = 1)
data("spectra", package = "MetCirc") ## create similarity matrix similarityMat <- Spectra::compareSpectra(sps_tissue[1:10], FUN = MsCoreUtils::ndotproduct, ppm = 20, m = 0.5, n = 2) rownames(similarityMat) <- colnames(similarityMat) <- sps_tissue$name[1:10] linkDf <- createLinkDf(similarityMatrix = similarityMat, sps = sps_tissue[1:10], condition = c("SPL", "LIM", "ANT", "STY"), lower = 0.01, upper = 1) ## cut link data.frame (here: only display links between groups) linkDf_cut <- cutLinkDf(linkDf, type = "inter") groupname <- c(as.character(linkDf_cut[, "spectrum1"]), as.character(linkDf_cut[, "spectrum2"])) groupname <- unique(groupname) ## plot legend circosLegend(groupname, highlight = TRUE, colour = NULL, cex = 1)
The 'data.frame' 'compartmentTissue' is used in the subsection 'Preparing the tissue data set for analysis' in the vignette of 'MetCirc'. In 'compartmentTissue', information on the organ-localisation of each MS/MS spectrum is stored.
'data.frame'
'data.frame'
Thomas Naake, [email protected]
internal
'convertExampleDF' is a 'data.frame' which comprises information on a specific metabolite per row stating the average retention time, average m/z, the name of the metabolite, the adduct ion name and the spectrum reference file name. The function 'allocatePrecursor2mz' uses 'data.frame's of the kind of 'sd01\_outputXCMS' and 'sd02\_deconvoluted' to create a 'data.frame' of the kind of 'convertExampleDF'. Allocation of precursor ions to candidate m/z values is based on minimal distance of m/z and deviance of retention time based on an objective function. See '?allocatePrecursor2mz' for further information.
'data.frame'
'data.frame'
Thomas Naake, [email protected]
internal
Spectra
Convert msp data frame into object of class [Spectra()]
convertMsp2Spectra(msp)
convertMsp2Spectra(msp)
msp |
|
msp is a data frame of a .msp file, a typical data file for MS/MS libraries. The data frame has two columns and contains in the first column the entries "NAME:", "PRECURSORMZ:" (or "EXACTMASS:"), "RETENTIONTIME:", Num Peaks:" and information on fragments and peak areas/intensities and will extract the respective information in the second column.
convertMsp2Spectra
returns an object of class 'Spectra'
Thomas Naake, [email protected]
data("convertMsp2Spectra", package = "MetCirc") convertMsp2Spectra(msp = msp2spectra)
data("convertMsp2Spectra", package = "MetCirc") convertMsp2Spectra(msp = msp2spectra)
Create a link matrix which links every feature in similarity matrix with another.
createLink0df(similarityMatrix, sps, condition)
createLink0df(similarityMatrix, sps, condition)
similarityMatrix |
|
sps |
|
condition |
|
createLink0df
creates a matrix
from a similarity matrix which
includes all connections between features in the similarity matrix, but
exclude links which have a similarity of exactly 0.
createLink0df
returns a 'matrix' that gives per each row
information on linked features
Thomas Naake, [email protected]
data("spectra", package = "MetCirc") data("similarityMat", package = "MetCirc") link0df <- createLink0df(similarityMatrix = similarityMat, sps = sps_tissue, condition = c("SPL", "LIM", "ANT", "STY"))
data("spectra", package = "MetCirc") data("similarityMat", package = "MetCirc") link0df <- createLink0df(similarityMatrix = similarityMat, sps = sps_tissue, condition = c("SPL", "LIM", "ANT", "STY"))
Create a data frame which contains features to link (indices).
createLinkDf(similarityMatrix, sps, condition, lower, upper)
createLinkDf(similarityMatrix, sps, condition, lower, upper)
similarityMatrix |
|
sps |
|
condition |
|
lower |
|
upper |
|
lower
and upper
are numerical values and truncate
spectra based on their similarity. The function createLinkDf
is
a wrapper for the functions createLink0df
and thresholdLinkDf
.
createLinkDf
returns a data.frame
that gives per each
row information on linked features
Thomas Naake, [email protected]
data("spectra", package = "MetCirc") data("similarityMat", package = "MetCirc") link0df <- createLink0df(similarityMatrix = similarityMat, sps = sps_tissue, condition = c("SPL", "LIM", "ANT", "STY")) createLinkDf(similarityMatrix = similarityMat, sps = sps_tissue, condition = c("SPL", "LIM", "ANT", "STY"), lower = 0.5, upper = 1)
data("spectra", package = "MetCirc") data("similarityMat", package = "MetCirc") link0df <- createLink0df(similarityMatrix = similarityMat, sps = sps_tissue, condition = c("SPL", "LIM", "ANT", "STY")) createLinkDf(similarityMatrix = similarityMat, sps = sps_tissue, condition = c("SPL", "LIM", "ANT", "STY"), lower = 0.5, upper = 1)
Create a cut link data frame
cutLinkDf(linkDf, type = c("all", "inter", "intra"))
cutLinkDf(linkDf, type = c("all", "inter", "intra"))
linkDf |
|
type |
|
This function is used to truncate features from linkDf
. If
type = "all"
, linkDf
will not be changed; if
type = "inter"
the returned linkDf
will only contain entries
of links which are between groups and not inside groups; contrary to that,
if type = "intra"
the returned linkDf
will only contain entries
of links which are inside groups and not between groups.
cutLinkDf
returns a data.frame
that gives per each row
information on linked features
Thomas Naake, [email protected]
data("spectra", package = "MetCirc") data("similarityMat", package = "MetCirc") linkDf <- createLinkDf(similarityMatrix = similarityMat, sps = sps_tissue, condition = c("SPL", "LIM", "ANT", "STY"), lower = 0.75, upper = 1) cutLinkDf(linkDf = linkDf, type = "all")
data("spectra", package = "MetCirc") data("similarityMat", package = "MetCirc") linkDf <- createLinkDf(similarityMatrix = similarityMat, sps = sps_tissue, condition = c("SPL", "LIM", "ANT", "STY"), lower = 0.75, upper = 1) cutLinkDf(linkDf = linkDf, type = "all")
linkDf
of featureGets indices in linkDf
of feature
getLinkDfIndices(groupnameselected, linkDf)
getLinkDfIndices(groupnameselected, linkDf)
groupnameselected |
|
linkDf |
|
Internal use for function highlight
getLinkDfIndices
returns indices concerning linkDf
to which
groupnameselected
connects
Thomas Naake, [email protected]
## Not run: getLinkDfIndices(groupnameselected, linkMatrix)
## Not run: getLinkDfIndices(groupnameselected, linkMatrix)
A function to add links and highlight sectors to an initialised
and plotted circlize
plot with one track.
highlight( groupname, ind, linkDf, colour = NULL, transparency = 0.4, links = TRUE )
highlight( groupname, ind, linkDf, colour = NULL, transparency = 0.4, links = TRUE )
groupname |
|
ind |
|
linkDf |
|
colour |
|
transparency |
|
links |
|
Internal use for shinyCircos
or outside of shinyCircos
to
reproduce the figure.
The function will update an existing plot by highlighting a specified sector and connected links.
Thomas Naake, [email protected]
data("spectra", package = "MetCirc") ## create similarity matrix similarityMat <- Spectra::compareSpectra(sps_tissue[1:10], FUN = MsCoreUtils::ndotproduct, ppm = 20, m = 0.5, n = 2) rownames(similarityMat) <- colnames(similarityMat) <- sps_tissue$name[1:10] ## order similarityMat according to retentionTime and update rownames simM <- orderSimilarityMatrix(similarityMat, sps = sps_tissue[1:10], type = "retentionTime") ## create link matrix linkDf <- createLinkDf(similarityMatrix = simM, sps = sps_tissue, condition = c("SPL", "LIM", "ANT", "STY"), lower = 0.01, upper = 1) ## cut link matrix (here: only display links between groups) linkDf_cut <- cutLinkDf(linkDf, type = "inter") ## set circlize parameters circos.clear() circos.par(gap.degree = 0, cell.padding = c(0.0, 0, 0.0, 0), track.margin = c(0.0, 0)) groupname <- c(as.character(linkDf_cut[, "spectrum1"]), as.character(linkDf_cut[, "spectrum2"])) groupname <- unique(groupname) ## here: set indSelected arbitrarily indSelected <- c(2,3) ## actual plotting plotCircos(groupname, linkDf_cut, initialize = TRUE, featureNames = TRUE, cexFeatureNames = 0.2, groupSector = TRUE, groupName = FALSE, links = FALSE, highlight = TRUE) ## highlight highlight(groupname = groupname, ind = indSelected, linkDf = linkDf_cut, colour = NULL, transparency = 0.4, links = TRUE)
data("spectra", package = "MetCirc") ## create similarity matrix similarityMat <- Spectra::compareSpectra(sps_tissue[1:10], FUN = MsCoreUtils::ndotproduct, ppm = 20, m = 0.5, n = 2) rownames(similarityMat) <- colnames(similarityMat) <- sps_tissue$name[1:10] ## order similarityMat according to retentionTime and update rownames simM <- orderSimilarityMatrix(similarityMat, sps = sps_tissue[1:10], type = "retentionTime") ## create link matrix linkDf <- createLinkDf(similarityMatrix = simM, sps = sps_tissue, condition = c("SPL", "LIM", "ANT", "STY"), lower = 0.01, upper = 1) ## cut link matrix (here: only display links between groups) linkDf_cut <- cutLinkDf(linkDf, type = "inter") ## set circlize parameters circos.clear() circos.par(gap.degree = 0, cell.padding = c(0.0, 0, 0.0, 0), track.margin = c(0.0, 0)) groupname <- c(as.character(linkDf_cut[, "spectrum1"]), as.character(linkDf_cut[, "spectrum2"])) groupname <- unique(groupname) ## here: set indSelected arbitrarily indSelected <- c(2,3) ## actual plotting plotCircos(groupname, linkDf_cut, initialize = TRUE, featureNames = TRUE, cexFeatureNames = 0.2, groupSector = TRUE, groupName = FALSE, links = FALSE, highlight = TRUE) ## highlight highlight(groupname = groupname, ind = indSelected, linkDf = linkDf_cut, colour = NULL, transparency = 0.4, links = TRUE)
Calculates the nearest feature in polar coordinates given cartesian coordinates.
minFragCart2Polar(x, y, degreeOfFeatures)
minFragCart2Polar(x, y, degreeOfFeatures)
x |
|
y |
|
degreeOfFeatures |
|
minFragCart2Polar
is employed to find the feature with
the smallest distance from given cartesian coordinates.
minFragCart2Polar
returns the index of the feature that has the
smallest distance to the given coordinates. As minFragCart2Polar
is
used in shinyCircos
for the track 1 only polar r
coordinates
between 0.8 and 1 will be used to find the feature with smallest distance.
Thomas Naake, [email protected]
library("MsCoreUtils") data("spectra", package = "MetCirc") ## create similarity matrix similarityMat <- Spectra::compareSpectra(sps_tissue[1:10], FUN = MsCoreUtils::ndotproduct, ppm = 20, m = 0.5, n = 2) rownames(similarityMat) <- colnames(similarityMat) <- sps_tissue$name[1:10] linkDf <- createLinkDf(similarityMatrix = similarityMat, sps = sps_tissue[1:10], condition = c("SPL", "LIM", "ANT", "STY"), lower = 0.5, upper = 1) ## cut link data.frame (here: only display links between groups) linkDf_cut <- cutLinkDf(linkDf, type = "inter") groupname <- c(as.character(linkDf_cut[, "spectrum1"]), as.character(linkDf_cut[, "spectrum2"])) groupname <- unique(groupname) ## set circlize parameters circos.clear() circos.par(gap.degree = 0, cell.padding = c(0.0, 0, 0.0, 0), track.margin = c(0.0, 0)) plotCircos(groupname, NULL, initialize = TRUE, featureNames = FALSE, groupName = FALSE, groupSector = FALSE, links = FALSE, highlight = FALSE) x <- 1 y <- 0 degreeFeatures <- lapply(groupname, function(x) mean(circlize:::get.sector.data(x)[c("start.degree", "end.degree")])) minFragCart2Polar(x, y, degreeOfFeatures = degreeFeatures)
library("MsCoreUtils") data("spectra", package = "MetCirc") ## create similarity matrix similarityMat <- Spectra::compareSpectra(sps_tissue[1:10], FUN = MsCoreUtils::ndotproduct, ppm = 20, m = 0.5, n = 2) rownames(similarityMat) <- colnames(similarityMat) <- sps_tissue$name[1:10] linkDf <- createLinkDf(similarityMatrix = similarityMat, sps = sps_tissue[1:10], condition = c("SPL", "LIM", "ANT", "STY"), lower = 0.5, upper = 1) ## cut link data.frame (here: only display links between groups) linkDf_cut <- cutLinkDf(linkDf, type = "inter") groupname <- c(as.character(linkDf_cut[, "spectrum1"]), as.character(linkDf_cut[, "spectrum2"])) groupname <- unique(groupname) ## set circlize parameters circos.clear() circos.par(gap.degree = 0, cell.padding = c(0.0, 0, 0.0, 0), track.margin = c(0.0, 0)) plotCircos(groupname, NULL, initialize = TRUE, featureNames = FALSE, groupName = FALSE, groupSector = FALSE, links = FALSE, highlight = FALSE) x <- 1 y <- 0 degreeFeatures <- lapply(groupname, function(x) mean(circlize:::get.sector.data(x)[c("start.degree", "end.degree")])) minFragCart2Polar(x, y, degreeOfFeatures = degreeFeatures)
'convertMsp2Spectra' contains the object 'msp2spectra' that is a data frame in .MSP format, a typical format for MS/MS library building. Each entry consists of the metabolite name (NAME), the precursor mz (PRECURSORMZ), the retention time (RETENTIONTIME), number of peaks (Num Peaks), together with fragments and their intensity values. In the example used in the function 'convertMsp2Spectra' the 'data.frame' 'msp2spectra' is used to construct an object of class 'MSpectra'.
'data.frame'
'data.frame'
Thomas Naake, [email protected]
http://prime.psc.riken.jp/Metabolomics_Software/MS-DIAL/, truncated .MSP file of GNPS MS/MS Negative (contains 22 entries): http://prime.psc.riken.jp/Metabolomics_Software/MS-DIAL/MSMS-GNPS-Curated-Neg.msp
Calculate similarity based on neutral losses (NLS)
neutralloss(x, y, m = 0.5, n = 2, na.rm = TRUE, ...)
neutralloss(x, y, m = 0.5, n = 2, na.rm = TRUE, ...)
x |
|
y |
|
m |
|
n |
|
na.rm |
|
... |
further arguments |
Similarities of spectra based on neutral losses are calculated according to the following formula:
,
with and
. For further information
see Li et al. (2015): Navigating natural variation in herbivory-induced
secondary metabolism in coyote tobacco populations using MS/MS structural
analysis. PNAS, E4147–E4155.
In here, the precursor m/z is taken by the m/z feature with the highest intensity.
neutralloss
returns a numeric value ranging between 0 and 1, where 0
indicates no similarity between the two MS/MS features, while 1 indicates
that the MS/MS features are identical.
Prior to calculating
or
, all intensity values are divided by the maximum intensity value.
neutralloss
returns a numeric similarity coefficient between 0 and 1
Thomas Naake, [email protected]
data("spectra", package = "MetCirc") Spectra::compareSpectra(sps_tissue[1:10], FUN = neutralloss, m = 0.5, n = 2)
data("spectra", package = "MetCirc") Spectra::compareSpectra(sps_tissue[1:10], FUN = neutralloss, m = 0.5, n = 2)
Internal function for shiny application. May also be used outside of shiny to reconstruct figures.
orderSimilarityMatrix( similarityMatrix, sps, type = c("retentionTime", "mz", "clustering"), group = FALSE )
orderSimilarityMatrix( similarityMatrix, sps, type = c("retentionTime", "mz", "clustering"), group = FALSE )
similarityMatrix |
|
sps |
|
type |
|
group |
|
orderSimilarityMatrix
takes a similarity matrix,
Spectra
object (sps
, containing information on m/z and
retention time), and a character
vector as arguments. It will then
reorder rows and columns of the similarityMatrix
object such,
that it orders rows and columns of similarityMatrix
according to
m/z, retention time or clustering in each group.
orderSimilarityMatrix
is employed in the shinyCircos
function to create similarityMatrix
objects which will allow to switch
between different types of ordering in between groups (sectors) in the
circos plot. It may be used as well externally, to reproduce plots outside
of the reactive environment (see vignette for a workflow).
matrix
, orderSimilarityMatrix
returns a similarity
matrix with ordered rownames
according to the character
vector
type
Thomas Naake, [email protected]
data("spectra", package = "MetCirc") similarityMat <- Spectra::compareSpectra(sps_tissue[1:10], FUN = MsCoreUtils::ndotproduct, ppm = 10) rownames(similarityMat) <- colnames(similarityMat) <- sps_tissue$name[1:10] ## order according to retention time orderSimilarityMatrix(similarityMatrix = similarityMat, sps = sps_tissue, type = "retentionTime", group = FALSE)
data("spectra", package = "MetCirc") similarityMat <- Spectra::compareSpectra(sps_tissue[1:10], FUN = MsCoreUtils::ndotproduct, ppm = 10) rownames(similarityMat) <- colnames(similarityMat) <- sps_tissue$name[1:10] ## order according to retention time orderSimilarityMatrix(similarityMatrix = similarityMat, sps = sps_tissue, type = "retentionTime", group = FALSE)
Circular plot to visualize similarity.
plotCircos( groupname, linkDf, initialize = c(TRUE, FALSE), featureNames = c(TRUE, FALSE), cexFeatureNames = 0.3, groupSector = c(TRUE, FALSE), groupName = c(TRUE, FALSE), links = c(TRUE, FALSE), highlight = c(TRUE, FALSE), colour = NULL, transparency = 0.2 )
plotCircos( groupname, linkDf, initialize = c(TRUE, FALSE), featureNames = c(TRUE, FALSE), cexFeatureNames = 0.3, groupSector = c(TRUE, FALSE), groupName = c(TRUE, FALSE), links = c(TRUE, FALSE), highlight = c(TRUE, FALSE), colour = NULL, transparency = 0.2 )
groupname |
|
linkDf |
|
initialize |
|
featureNames |
|
cexFeatureNames |
|
groupSector |
|
groupName |
|
links |
|
highlight |
|
colour |
|
transparency |
|
Internal use for shinyCircos
or used outside of
shinyCircos
to reproduce figure
The function will initialize a circlize plot and/or will plot features of a circlize plot.
Thomas Naake, [email protected]
library("MsCoreUtils") data("spectra", package = "MetCirc") ## create similarity matrix similarityMat <- Spectra::compareSpectra(sps_tissue[1:10], FUN = MsCoreUtils::ndotproduct, ppm = 20, m = 0.5, n = 2) rownames(similarityMat) <- colnames(similarityMat) <- sps_tissue$name[1:10] ## order similarityMat according to retentionTime simM <- orderSimilarityMatrix(similarityMat, sps = sps_tissue[1:10], type = "retentionTime") ## create link data.frame linkDf <- createLinkDf(similarityMatrix = simM, sps = sps_tissue, condition = c("SPL", "LIM", "ANT", "STY"), lower = 0.01, upper = 1) ## cut link data.frame (here: only display links between groups) linkDf_cut <- cutLinkDf(linkDf, type = "inter") ## set circlize paramters circos.clear() circos.par(gap.degree = 0, cell.padding = c(0.0, 0, 0.0, 0), track.margin = c(0.0, 0)) groupname <- c(as.character(linkDf_cut[, "spectrum1"]), as.character(linkDf_cut[, "spectrum2"])) groupname <- unique(groupname) ## actual plotting plotCircos(groupname, linkDf_cut, initialize = TRUE, featureNames = TRUE, cexFeatureNames = 0.3, groupSector = TRUE, groupName = FALSE, links = FALSE, highlight = FALSE, colour = NULL, transparency = 0.2)
library("MsCoreUtils") data("spectra", package = "MetCirc") ## create similarity matrix similarityMat <- Spectra::compareSpectra(sps_tissue[1:10], FUN = MsCoreUtils::ndotproduct, ppm = 20, m = 0.5, n = 2) rownames(similarityMat) <- colnames(similarityMat) <- sps_tissue$name[1:10] ## order similarityMat according to retentionTime simM <- orderSimilarityMatrix(similarityMat, sps = sps_tissue[1:10], type = "retentionTime") ## create link data.frame linkDf <- createLinkDf(similarityMatrix = simM, sps = sps_tissue, condition = c("SPL", "LIM", "ANT", "STY"), lower = 0.01, upper = 1) ## cut link data.frame (here: only display links between groups) linkDf_cut <- cutLinkDf(linkDf, type = "inter") ## set circlize paramters circos.clear() circos.par(gap.degree = 0, cell.padding = c(0.0, 0, 0.0, 0), track.margin = c(0.0, 0)) groupname <- c(as.character(linkDf_cut[, "spectrum1"]), as.character(linkDf_cut[, "spectrum2"])) groupname <- unique(groupname) ## actual plotting plotCircos(groupname, linkDf_cut, initialize = TRUE, featureNames = TRUE, cexFeatureNames = 0.3, groupSector = TRUE, groupName = FALSE, links = FALSE, highlight = FALSE, colour = NULL, transparency = 0.2)
plotSpectra
plots a spectra of a subject
and
query
spectra. plotSpectra
uses ggplot
plotting functionality.
plotSpectra(sps, subject, query)
plotSpectra(sps, subject, query)
sps |
|
subject |
|
query |
|
Internally, all intensities are normalized to 100%.
ggplot2
plot
Thomas Naake, [email protected]
data("spectra", package = "MetCirc") plotSpectra(sps = sps_tissue, subject = "SPL_1", query = "SPL_2")
data("spectra", package = "MetCirc") plotSpectra(sps = sps_tissue, subject = "SPL_1", query = "SPL_2")
Displays information on connected features of selected features.
printInformationSelect( select, sps = NULL, linkDfInd, linkDf, similarityMatrix, roundDigits = 2 )
printInformationSelect( select, sps = NULL, linkDfInd, linkDf, similarityMatrix, roundDigits = 2 )
select |
|
sps |
|
linkDfInd |
|
linkDf |
|
similarityMatrix |
|
roundDigits |
|
printInformationSelect
is for internal use.
character
that is in HTML format
Thomas Naake, [email protected]
data("spectra", package = "MetCirc") sps_tissue@metadata$names <- rep("Unknown", 259) sps_tissue@metadata$information <- rep("Unknown", 259) sps_tissue@metadata$classes <- rep("Unknown", 259) sps_tissue@metadata$adduct <- rep("Unknown", 259) similarityMat <- Spectra::compareSpectra(sps_tissue[1:10], FUN = MsCoreUtils::ndotproduct, ppm = 20, m = 0.5, n = 2) rownames(similarityMat) <- colnames(similarityMat) <- sps_tissue$name[1:10] linkDf <- createLinkDf(similarityMatrix = similarityMat, sps = sps_tissue[1:10], condition = c("SPL", "LIM", "ANT", "STY"), lower = 0.01, upper = 1) ## cut link data.frame (here: only display links between groups) linkDf_cut <- cutLinkDf(linkDf, type = "inter") groupname <- c(as.character(linkDf_cut[, "spectrum1"]), as.character(linkDf_cut[, "spectrum2"])) groupname <- unique(groupname) ## arbitrarily select a feature ind <- 2 linkDfInds <- getLinkDfIndices(groupname[ind], linkDf_cut) MetCirc:::printInformationSelect(select = groupname[ind], sps = sps_tissue[1:10], linkDfInd = linkDfInds, linkDf = linkDf_cut, similarityMatrix = similarityMat)
data("spectra", package = "MetCirc") sps_tissue@metadata$names <- rep("Unknown", 259) sps_tissue@metadata$information <- rep("Unknown", 259) sps_tissue@metadata$classes <- rep("Unknown", 259) sps_tissue@metadata$adduct <- rep("Unknown", 259) similarityMat <- Spectra::compareSpectra(sps_tissue[1:10], FUN = MsCoreUtils::ndotproduct, ppm = 20, m = 0.5, n = 2) rownames(similarityMat) <- colnames(similarityMat) <- sps_tissue$name[1:10] linkDf <- createLinkDf(similarityMatrix = similarityMat, sps = sps_tissue[1:10], condition = c("SPL", "LIM", "ANT", "STY"), lower = 0.01, upper = 1) ## cut link data.frame (here: only display links between groups) linkDf_cut <- cutLinkDf(linkDf, type = "inter") groupname <- c(as.character(linkDf_cut[, "spectrum1"]), as.character(linkDf_cut[, "spectrum2"])) groupname <- unique(groupname) ## arbitrarily select a feature ind <- 2 linkDfInds <- getLinkDfIndices(groupname[ind], linkDf_cut) MetCirc:::printInformationSelect(select = groupname[ind], sps = sps_tissue[1:10], linkDfInd = linkDfInds, linkDf = linkDf_cut, similarityMatrix = similarityMat)
recordPlotFill_degreeFeatures
records a plot of filled
features and returns the degree of features.
recordPlotFill_degreeFeatures(type_match, ...)
recordPlotFill_degreeFeatures(type_match, ...)
type_match |
|
... |
further arguments passed to |
Helper function for shinyCircos
.
list
of length 2, entry plotFill
is of recordedplot
and
entry degreeFeatures
is a list
of vectors of numeric(1)
Thomas Naake, [email protected]
type_match <- c("a_1", "a_2", "a_3", "b_1", "b_2", "b_3", "c_1", "c_2") MetCirc:::recordPlotFill_degreeFeatures(type_match)
type_match <- c("a_1", "a_2", "a_3", "b_1", "b_2", "b_3", "c_1", "c_2") MetCirc:::recordPlotFill_degreeFeatures(type_match)
recordedplot
of plotCircos
plot with
highlight = TRUE
recordPlotHighlight
returns a recordedplot
object of
plotCircos
with highlight = TRUE
recordPlotHighlight(type_match, ...)
recordPlotHighlight(type_match, ...)
type_match |
|
... |
further arguments passed to |
Helper function for shinyCircos
.
recordedplot
Thomas Naake, [email protected]
type_match <- c("a_1", "a_2", "a_3", "b_1", "b_2", "b_3", "c_1", "c_2") MetCirc:::recordPlotHighlight(type_match)
type_match <- c("a_1", "a_2", "a_3", "b_1", "b_2", "b_3", "c_1", "c_2") MetCirc:::recordPlotHighlight(type_match)
replayPlotAdd
plots additional plots on a plot, either
plots plotCircos
or highlight
.
replayPlotAdd( orderMatch = "mz", onCircle = FALSE, linkDf, mz_match, rt_match, clust_match, ind, indMz, indRT, indCluster )
replayPlotAdd( orderMatch = "mz", onCircle = FALSE, linkDf, mz_match, rt_match, clust_match, ind, indMz, indRT, indCluster )
orderMatch |
|
onCircle |
|
linkDf |
|
mz_match |
|
rt_match |
|
clust_match |
|
ind |
|
indMz |
|
indRT |
|
indCluster |
|
Helper function for shinyCircos
.
Depending on onCircle
and indMz
either returns
plotCircos
or highlight
Thomas Naake, [email protected]
data("spectra", package = "MetCirc") similarityMat <- Spectra::compareSpectra(sps_tissue[1:10], FUN = MsCoreUtils::ndotproduct, ppm = 10, m = 0.5, n = 2) rownames(similarityMat) <- colnames(similarityMat) <- sps_tissue$name[1:10] ## order according to m/z mz_match <- MetCirc:::typeMatch_link0(similarityMatrix = similarityMat, sps = sps_tissue, type = "mz", condition = c("SPL", "LIM", "ANT", "STY")) linkDf <- mz_match[["link0df"]] mz_match <- mz_match[["type_match"]] ## order according to retention time rt_match <- MetCirc:::typeMatch_link0(similarityMatrix = similarityMat, sps = sps_tissue, type = "retentionTime", condition = c("SPL", "LIM", "ANT", "STY")) rt_match <- rt_match[["type_match"]] ## order according to clustering clust_match <- MetCirc:::typeMatch_link0(similarityMatrix = similarityMat, sps = sps_tissue, type = "clustering", condition = c("SPL", "LIM", "ANT", "STY")) clust_match <- clust_match[["type_match"]] circos.initialize(mz_match,##, levels = mz_match), xlim = matrix(rep(c(0,1), length(mz_match)), ncol = 2, byrow = TRUE)) #circos.trackPlotRegion(factor(mz_match, levels = mz_match), ylim = c(0,1)) MetCirc:::replayPlotAdd(orderMatch = "mz", onCircle = FALSE, linkDf = linkDf, mz_match = mz_match, rt_match = rt_match, clust_match = clust_match, ind = 1, indMz = NULL, indRT = NULL, indCluster = NULL)
data("spectra", package = "MetCirc") similarityMat <- Spectra::compareSpectra(sps_tissue[1:10], FUN = MsCoreUtils::ndotproduct, ppm = 10, m = 0.5, n = 2) rownames(similarityMat) <- colnames(similarityMat) <- sps_tissue$name[1:10] ## order according to m/z mz_match <- MetCirc:::typeMatch_link0(similarityMatrix = similarityMat, sps = sps_tissue, type = "mz", condition = c("SPL", "LIM", "ANT", "STY")) linkDf <- mz_match[["link0df"]] mz_match <- mz_match[["type_match"]] ## order according to retention time rt_match <- MetCirc:::typeMatch_link0(similarityMatrix = similarityMat, sps = sps_tissue, type = "retentionTime", condition = c("SPL", "LIM", "ANT", "STY")) rt_match <- rt_match[["type_match"]] ## order according to clustering clust_match <- MetCirc:::typeMatch_link0(similarityMatrix = similarityMat, sps = sps_tissue, type = "clustering", condition = c("SPL", "LIM", "ANT", "STY")) clust_match <- clust_match[["type_match"]] circos.initialize(mz_match,##, levels = mz_match), xlim = matrix(rep(c(0,1), length(mz_match)), ncol = 2, byrow = TRUE)) #circos.trackPlotRegion(factor(mz_match, levels = mz_match), ylim = c(0,1)) MetCirc:::replayPlotAdd(orderMatch = "mz", onCircle = FALSE, linkDf = linkDf, mz_match = mz_match, rt_match = rt_match, clust_match = clust_match, ind = 1, indMz = NULL, indRT = NULL, indCluster = NULL)
replayPlot
replayPlotOrder
will call replayPlot
from grDevices
with
a recordedplot
object based on orderMatch
.
replayPlotOrder(orderMatch = "mz", onCircle = FALSE, plot_l, ind)
replayPlotOrder(orderMatch = "mz", onCircle = FALSE, plot_l, ind)
orderMatch |
|
onCircle |
|
plot_l |
|
ind |
|
Helper function for shinyCircos
.
replayedplot
Thomas Naake, [email protected]
type_match <- c("a_1", "a_2", "a_3", "b_1", "b_2", "b_3", "c_1", "c_2") plotCircos(type_match, NULL, initialize = TRUE, featureNames = TRUE, groupSector = TRUE, groupName = FALSE, links = FALSE, highlight = TRUE) p <- recordPlot() plot.new() plot_l <- list(highlightMz = p) MetCirc:::replayPlotOrder(orderMatch = "mz", onCircle = TRUE, plot_l = plot_l, ind = NULL)
type_match <- c("a_1", "a_2", "a_3", "b_1", "b_2", "b_3", "c_1", "c_2") plotCircos(type_match, NULL, initialize = TRUE, featureNames = TRUE, groupSector = TRUE, groupName = FALSE, links = FALSE, highlight = TRUE) p <- recordPlot() plot.new() plot_l <- list(highlightMz = p) MetCirc:::replayPlotOrder(orderMatch = "mz", onCircle = TRUE, plot_l = plot_l, ind = NULL)
'sd01_outputXCMS' is the output file from the package 'XCMS' using the data from Li et al. (2015). See Li et al. (2015) for further details.
'data.frame'
'data.frame'
Thomas Naake, [email protected]
Li et al. (2015)
'sd02_deconvoluted' contains MS/MS data from Li et al. (2015). It is a 'data.frame' which hosts m/z values, retention time, intensity and the respective precursor m/z values. 'sd02_deconvoluted' originates from Li et al. (2015). See Li et al. (2015) for further information.
'data.frame'
'data.frame'
Thomas Naake, [email protected]
Li et al. (2015)
select
returns mz
, rt
or clust
depending on
condition
.
select(condition, mz, rt, clust)
select(condition, mz, rt, clust)
condition |
|
mz |
object to return if |
rt |
object to return if |
clust |
object to return if |
Helper function for shinyCircos
, replayPlotOrder
and
replayPlotAdd
.
mz
, rt
or clust
depending on condition
Thomas Naake [email protected]
mz <- 1 rt <- 2 clust <- 3 MetCirc:::select(condition = "mz", mz = mz, rt = rt, clust = clust)
mz <- 1 rt <- 2 clust <- 3 MetCirc:::select(condition = "mz", mz = mz, rt = rt, clust = clust)
Visualise the similarity of MS/MS features in a reactive
context. See Details
the vignette for further descriptions on how to
use shinyCircos
.
shinyCircos(similarityMatrix, sps, condition, ...)
shinyCircos(similarityMatrix, sps, condition, ...)
similarityMatrix |
|
sps |
|
condition |
|
... |
further arguments passed to |
The function is based on the shiny
and circlize
package.
The user can choose interactively thresholds, type of links (between or
within groups), display information about MS/MS features, permanently select
MS/MS features and export selected precursors. The Spectra
object
stores annotation information about the MS/MS features. Names of features
within the similarityMatrix
have to be found as entries
in Spectra
. sps$name
are used as identifiers and
colnames
/rownames
from similarityMatrix
are cleaved
by the group identifier (separated by "_"
). Annotation information
is taken from spectra
from the columns "names", "information",
"classes" and "adduct" in the slot metadata
of spectra
.
After exiting the application, the annotation will be written to the
respective columns in the slot metadata
. If one or several of these
columns is already present in metadata
, the column(s) will be used as
the source of annotation information.
character
, shinyCircos
returns a character
vector with the permanently selected precursors and an object with the
MSpectra
object containing the annotation.
Thomas Naake, [email protected]
data("spectra", package = "MetCirc") similarityMat <- Spectra::compareSpectra(sps_tissue[1:10], FUN = MsCoreUtils::ndotproduct, ppm = 10, m = 0.5, n = 2) rownames(similarityMat) <- colnames(similarityMat) <- sps_tissue$name[1:10] ## Not run: shinyCircos(similarityMatrix = similarityMat, sps = sps_tissue, condition = c("SPL", "LIM", "ANT", "STY")) ## End(Not run)
data("spectra", package = "MetCirc") similarityMat <- Spectra::compareSpectra(sps_tissue[1:10], FUN = MsCoreUtils::ndotproduct, ppm = 10, m = 0.5, n = 2) rownames(similarityMat) <- colnames(similarityMat) <- sps_tissue$name[1:10] ## Not run: shinyCircos(similarityMatrix = similarityMat, sps = sps_tissue, condition = c("SPL", "LIM", "ANT", "STY")) ## End(Not run)
'similarityMat' is a 'matrix' containing the pair-wise similarity scores derived from the 'idMSMStissueproject' data set. See the vignette for a workflow to reproduce the object 'similarityMat'.
'matrix'
'matrix'
Thomas Naake, [email protected]
data("spectra", package = "MetCirc") similarityMat <- Spectra::compareSpectra(sps_tissue, fun = ndotproduct, ppm = 10) rownames(similarityMat) <- colnames(similarityMat) <- sps_tissue$name save(similarityMat, file = "similarityMat.RData", compress = "xz")
spectraCondition
returns the names of Spectra
that are
present in condition, corresponding to the slot
metadata
.
spectraCondition(sps, condition)
spectraCondition(sps, condition)
sps |
|
condition |
|
Helper function in createLink0df
and shinyCircos
.
list
, named list
with character
vector as entries that
contains the names of the MS/MS entries in spectra
that are present
in the conditon
(tissues, stress conditions, time points, etc.)
Thomas Naake [email protected]
data("spectra", package = "MetCirc") MetCirc:::spectraCondition(sps = sps_tissue, condition = c("SPL", "LIM", "ANT", "STY"))
data("spectra", package = "MetCirc") MetCirc:::spectraCondition(sps = sps_tissue, condition = c("SPL", "LIM", "ANT", "STY"))
'sps_tissue' is a 'Spectra' object derived from the 'idMSMStissueproject' data set. See the vignette for a workflow to reproduce the object 'spectra'.
'matrix'
'matrix'
Thomas Naake, [email protected]
data("idMSMStissueproject", package = "MetCirc") ## get all MS/MS spectra tissue <- tissue[tissue[, "id"] id_uniq <- unique(tissue[, "id"])
## obtain precursor m/z from id_uniq prec_mz <- lapply(strsplit(as.character(id_uniq), split = "_"), "[", 1) |> unlist() |> as.numeric()
## obtain m/z from fragments per precursor m/z mz_l <- lapply(id_uniq, function(id_i) tissue[tissue[, "id"] == id_i, "mz"])
## obtain corresponding intensity values int_l <- lapply(id_uniq, function(id_i) tissue[tissue[, "id"] == id_i, "intensity"])
## order mz and intensity values int_l <- lapply(seq_along(int_l), function(i) int_l[[i]][order(mz_l[[i]])]) mz_l <- lapply(seq_along(mz_l), function(i) sort(mz_l[[i]]))
## obtain retention time by averaging all retention time values rt <- lapply(id_uniq, function(id_i) tissue[tissue[, "id"] == id_i, "rt"]) |> lapply(mean) |> unlist()
## create list of Spectra objects and concatenate sps_l <- lapply(seq_len(length(mz_l)), function(i) spd <- S4Vectors::DataFrame( name = as.character(i), rtime = rt[i], msLevel = 2L, precursorMz = prec_mz[i]) spd$mz <- list(mz_l[[i]]) spd$intensity <- list(int_l[[i]]) Spectra::Spectra(spd)) sps_tissue <- Reduce(c, sps_l)
## combine list of spectrum2 objects to MSpectra object, ## use SPL, LIM, ANT, STY for further analysis sps_tissue@metadata <- data.frame( compartmentTissue[, c("SPL", "LIM", "ANT", "STY")])
save(sps_tissue, file = "spectra.RData", compress = "xz")
Threshold a link data frame based on lower and upper similarity values. The function will return the links that lie within the defined bounds.
thresholdLinkDf(link0df, lower = 0.75, upper = 1)
thresholdLinkDf(link0df, lower = 0.75, upper = 1)
link0df |
|
lower |
|
upper |
|
lower
and upper
are numerical values
and truncate mass spectra based on their similarity values.
thresholdLinkDf
returns a data.frame
that gives per each row
information on linked features which are linked within certain thresholds.
Thomas Naake, [email protected]
data("spectra", package = "MetCirc") data("similarityMat", package = "MetCirc") link0df <- createLink0df(similarityMatrix = similarityMat, sps_tissue, condition = c("SPL", "LIM", "ANT", "STY")) thresholdLinkDf(link0df = link0df, lower = 0.5, upper = 1)
data("spectra", package = "MetCirc") data("similarityMat", package = "MetCirc") link0df <- createLink0df(similarityMatrix = similarityMat, sps_tissue, condition = c("SPL", "LIM", "ANT", "STY")) thresholdLinkDf(link0df = link0df, lower = 0.5, upper = 1)
The 'data.frame' 'tissue' is used in the subsection 'Preparing the tissue data set for analysis' in the vignette of 'MetCirc'. MS/MS data are merged across floral organs in this 'data.frame'.
'data.frame'
'data.frame'
Thomas Naake, [email protected]
internal
typeMatch_link0
returns a list with accessors "link0df"
and
"type_match"
typeMatch_link0(similarityMatrix, sps, type, condition)
typeMatch_link0(similarityMatrix, sps, type, condition)
similarityMatrix |
|
sps |
|
type |
|
condition |
|
Helper function for shinyCircos
.
list
of length 2, entry link0df
is a data.frame
and
entry type_match
is a character
vector
Thomas Naake, [email protected]
data("spectra", package = "MetCirc") similarityMat <- Spectra::compareSpectra(sps_tissue[1:10], FUN = MsCoreUtils::ndotproduct, ppm = 10, m = 0.5, n = 2) rownames(similarityMat) <- colnames(similarityMat) <- sps_tissue$name[1:10] ## order according to retention time MetCirc:::typeMatch_link0(similarityMatrix = similarityMat, sps = sps_tissue, type = "mz", condition = c("SPL", "LIM", "ANT", "STY"))
data("spectra", package = "MetCirc") similarityMat <- Spectra::compareSpectra(sps_tissue[1:10], FUN = MsCoreUtils::ndotproduct, ppm = 10, m = 0.5, n = 2) rownames(similarityMat) <- colnames(similarityMat) <- sps_tissue$name[1:10] ## order according to retention time MetCirc:::typeMatch_link0(similarityMatrix = similarityMat, sps = sps_tissue, type = "mz", condition = c("SPL", "LIM", "ANT", "STY"))