metapone | R Documentation |
Metapone conducts pathway tests for untargeted metabolomics data. It has three main characteristics: (1) expanded database combining SMPDB and Mummichog databases, with manual cleaning to remove redundancies; (2) A new weighted testing scheme to address the issue of metabolite-feature matching uncertainties; (3) Can consider positive mode and negative mode data in a single analysis.
metapone(dat=NULL, type=NULL, pa, hmdbCompMZ, pos.adductlist = c("M+H", "M+NH4", "M+Na", "M+ACN+H","M+ACN+Na", "M+2ACN+H", "2M+H", "2M+Na", "2M+ACN+H"), neg.adductlist = c("M-H","M-2H","M-2H+Na","M-2H+K", "M-2H+NH4","M-H2O-H","M-H+Cl", "M+Cl", "M+2Cl"), use.fractional.count=TRUE, match.tol.ppm=5, p.threshold=0.05, n.permu=200, fractional.count.power=0.5, max.match.count=10, use.fgsea = FALSE, use.meta = FALSE)
dat |
The list of test results. An element in the list should be postive ion mode test results or negative ion mode test results with four columns: m/z, retention time, p-value, test statistic. The package doesn't require both pos and neg to be present. One ion mode result is sufficient. Multiple ion mode results are allowed. |
type |
The list of corresponding ion mode of each element in dat. Each element in the list should be "pos" or "neg". The size of type should be consistent with the size of dat. |
pa |
Pathway information. A data frame with five columns: database pathway ID, pathway name, HMDB ID, KEGG ID, category of pathway. |
hmdbCompMZ |
the m/z values of common adduct ions of HMDB metaboites. See the help file of hmdbCompMZ for details. |
pos.adductlist |
The vector of positive adduct ions to be considered. |
neg.adductlist |
The vector of negative adduct ions to be considered. |
use.fractional.count |
A lot of features match to multiple metabolites by m/z. Whether to discount such matches by using fractional counts. |
match.tol.ppm |
The ppm level when conducting m/z match. |
p.threshold |
The threshold of p-values of metabolic features to be considered significant. |
n.permu |
The number of permutations in permutation test. |
fractional.count.power |
The fractional counts are taken to this power to transform the weights. |
max.match.count |
When calculating fractional counts, some features might be matched to too many. In that case the number of matches is capped by the value of max.match.count. |
use.fgsea |
Whether to use a GSEA type test when performing pathway testing. When it is FALSE, a permutation-based weighted hypergeometric test is performed. |
use.meta |
Whether to perform a GSEA type test with weighted metabolites. When it is FALSE, a GSEA type test is performed on weighted features. |
The method returns a generic S4 object of class "metapone.result":
@test.results |
A matrix with 8 columns: "p_value", "n_significant metabolites", "n_mapped_metabolites", "n_metabolites", "significant metabolites", "mapped_metabolites", "lfdr", "adjust.p". Each row is for a pathway. When using GSEA test, "ES", "NES", "nMoreExtreme" are returned additionally. |
@mapped.features |
A list. Each item is for a pathway. The item lists matched significant metabolites. |
The columns in test.result are the following:
p_value |
The p-value for each enrichment. |
n_significant metabolites |
The number of weighted significant metabolites associated with the enrichment. |
n_mapped_metabolites |
The number of weighted metabolites associated with the enrichment. |
n_metabolites |
The number of metabolites associated with the enrichment. |
significant metabolites |
A string with the names of significant metabolites that drive the enrichment. |
mapped_metabolites |
A string with the names of metabolites that drive the enrichment. |
lfdr |
The local fdr value for each enrichment. |
adjust.p |
The enrichment BH-adjusted p-value for each enrichment. |
ES |
The enrichment score (Avaliable in GSEA test). |
NES |
The enrichment score normalized to mean enrichment of random samples of the same size (Avaliable in GSEA test). |
nMoreExtreme |
The number of times a random metabolite set had a more extreme enrichment score value (Avaliable in GSEA test). |
Tianwei Yu (yutianwei@cuhk.edu.cn) Leqi Tian (leqitian@link.cuhk.edu.cn)
Small Molecule Pathway Database
pa
, hmdbCompMZ
data(hmdbCompMZ.metapone) data(pa) data(pos) data(neg) dat <- list(pos, neg) type <- list("pos", "neg") # Permutation-based weighted hypergeometric test r<-metapone(dat, type, pa, hmdbCompMZ=hmdbCompMZ.metapone, p.threshold=0.05, n.permu=100,fractional.count.power=0.5, max.match.count=10) hist(ptable(r)[,1]) # Metabolites based GSEA test r<-metapone(dat, type, pa, hmdbCompMZ=hmdbCompMZ.metapone, p.threshold=0.05, n.permu=100,fractional.count.power=0.5, max.match.count=10, use.fgsea = TRUE, use.meta = TRUE) hist(ptable(r)[,1]) # Features based GSEA test r<-metapone(dat, type, pa, hmdbCompMZ=hmdbCompMZ.metapone, p.threshold=0.05, n.permu=100,fractional.count.power=0.5, max.match.count=10, use.fgsea = FALSE, use.meta = FALSE) hist(ptable(r)[,1])
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.