Title: | Finding an Active Metabolic Module in Atom Transition Network |
---|---|
Description: | This package implements a metabolic network analysis pipeline to identify an active metabolic module based on high throughput data. The pipeline takes as input transcriptional and/or metabolic data and finds a metabolic subnetwork (module) most regulated between the two conditions of interest. The package further provides functions for module post-processing, annotation and visualization. |
Authors: | Anastasiia Gainullina [aut], Mariia Emelianova [aut], Alexey Sergushichev [aut, cre] |
Maintainer: | Alexey Sergushichev <[email protected]> |
License: | MIT + file LICENCE |
Version: | 1.5.0 |
Built: | 2024-12-10 06:20:17 UTC |
Source: | https://github.com/bioc/gatom |
Abbreviate lipid labels for lipid module
abbreviateLabels(module, orig.names, abbrev.names)
abbreviateLabels(module, orig.names, abbrev.names)
module |
Module to prepare |
orig.names |
whether to use original names from the dataset |
abbrev.names |
whether to use abbreviated names for all lipids |
module object with abbreviated labels
Add reactions without highly changing genes but with high average expression
addHighlyExpressedEdges(m, g, top = 3000)
addHighlyExpressedEdges(m, g, top = 3000)
m |
Metabolic module |
g |
Scored graph |
top |
Maximum rank value for the gene to be considered highly expressed |
module with added edges that correspond to high average expression
data(mEx) data(gEx) m <- addHighlyExpressedEdges(m = mEx, g = gEx)
data(mEx) data(gEx) m <- addHighlyExpressedEdges(m = mEx, g = gEx)
Collapse atoms belonging to the same metabolite into one vertex
collapseAtomsIntoMetabolites(m)
collapseAtomsIntoMetabolites(m)
m |
Metabolic module |
module in which atoms of the same metabolite are collapsed into one
data(mEx) m <- collapseAtomsIntoMetabolites(m = mEx)
data(mEx) m <- collapseAtomsIntoMetabolites(m = mEx)
Connect atoms belonging to the same metabolite with edges
connectAtomsInsideMetabolite(m)
connectAtomsInsideMetabolite(m)
m |
Metabolic module |
module in which atoms of the same metabolite are connected
data(mEx) m <- connectAtomsInsideMetabolite(m = mEx)
data(mEx) m <- connectAtomsInsideMetabolite(m = mEx)
Creates shinyCyJS widget from module
createShinyCyJSWidget( module, layout = list(name = "cose-bilkent", animate = FALSE, randomize = FALSE, nodeDimensionsIncludeLabels = TRUE), ... )
createShinyCyJSWidget( module, layout = list(name = "cose-bilkent", animate = FALSE, randomize = FALSE, nodeDimensionsIncludeLabels = TRUE), ... )
module |
Module |
layout |
Layout for the module |
... |
Other parameters |
html widget of input module
data(mEx) hw <- createShinyCyJSWidget(module = mEx)
data(mEx) hw <- createShinyCyJSWidget(module = mEx)
This package implements a metabolic network analysis pipeline to identify an active metabolic module based on high throughput data. The pipeline takes as input transcriptional and/or metabolic data and finds a metabolic subnetwork (module) most regulated between the two conditions of interest. The package further provides functions for module post-processing, annotation and visualization.
Data preprocessing:
prepareDE
, getMetDEMeta
, getGeneDEMeta
Graph creation:
makeMetabolicGraph
Graph scoring:
scoreGraph
Module postprocessing:
collapseAtomsIntoMetabolites
, connectAtomsInsideMetabolite
,
addHighlyExpressedEdges
, abbreviateLabels
Plotting module:
createShinyCyJSWidget
Exporting module:
saveModuleToHtml
, saveModuleToDot
,
saveModuleToPdf
, saveModuleToXgmml
For detailed pipeline analysis, see gatom vignette:
vignette("gatom-tutorial", package = "gatom")
Example data provided by gatom consists of:
metabolite differential abundance data (met.de.rawEx
),
gene differential expression data (gene.de.rawEx
),
KEGG-based network object (networkEx
),
KEGG-based metabolite database object (met.kegg.dbEx
),
Example organism annotation object (org.Mm.eg.gatom.annoEx
),
metabolic graph with atom topology (gEx
),
scored metabolic graph with atom topology (gsEx
),
and metabolic module (mEx
).
See file https://github.com/ctlab/gatom/blob/master/inst/scripts/example.R for details.
tibble/data.frame object
Default values for all columns are NULL which mean they are determined automatically.
getGeneDEMeta( gene.de.raw, org.gatom.anno, idColumn = NULL, idType = NULL, pvalColumn = NULL, logPvalColumn = NULL, log2FCColumn = NULL, baseMeanColumn = NULL, signalColumn = NULL, signalRankColumn = NULL )
getGeneDEMeta( gene.de.raw, org.gatom.anno, idColumn = NULL, idType = NULL, pvalColumn = NULL, logPvalColumn = NULL, log2FCColumn = NULL, baseMeanColumn = NULL, signalColumn = NULL, signalRankColumn = NULL )
gene.de.raw |
A table with differential expression results, an object convertable to data.frame. |
org.gatom.anno |
Organsim-specific annotation obtained from makeOrgGatomAnnotation function. |
idColumn |
Specifies column name with gene identifiers. |
idType |
Specifies type of gene IDs (one of the supported by annotation). |
pvalColumn |
Specifies column with p-values. |
logPvalColumn |
Specifies column with log p-values, if there is no such column one will be generated automatically. |
log2FCColumn |
Specifies column with log2-fold changes. |
baseMeanColumn |
Specifies column with average expression across samples. |
signalColumn |
Specifies column with identifier of the measured entity (such as gene ID for RNA-seq and probe ID for microarrays). Could be NULL (automatic, set from based on pval and log2FC columns), character (column name), or function (evaluated in a scope of original data frame) |
signalRankColumn |
Specifies how the genes are ranked from highly to lowly expressed, used in 'addHighlyExpressedEdgues' function. Could be NULL (automatic), character (column name) function (evaluated in a scope of original data frame). |
object with prepared columns for the analysis for gene data
data("org.Mm.eg.gatom.annoEx") data("gene.de.rawEx") de.meta <- getGeneDEMeta(gene.de.rawEx, org.gatom.anno = org.Mm.eg.gatom.annoEx)
data("org.Mm.eg.gatom.annoEx") data("gene.de.rawEx") de.meta <- getGeneDEMeta(gene.de.rawEx, org.gatom.anno = org.Mm.eg.gatom.annoEx)
Generate list of metabolic pathways from Reactome and KEGG databases
getMetabolicPathways( universe, metGenes, keggOrgCode, threshold = 0.01, includeReactome = TRUE, includeKEGG = TRUE )
getMetabolicPathways( universe, metGenes, keggOrgCode, threshold = 0.01, includeReactome = TRUE, includeKEGG = TRUE )
universe |
list of genes |
metGenes |
list of metabolic genes |
keggOrgCode |
KEGG organism code, like mmu or hsa |
threshold |
threshold for Fisher test to filter out non-metabolic pathways |
includeReactome |
whether to include Reactome pathways (only works for Entrez ID universe) |
includeKEGG |
whether to include KEGG pathways and modules |
list of metabolic pathways for given organism and list of genes
Finds columns in differential expression table for metabolites required for gatom analysis
getMetDEMeta( met.de.raw, met.db, idColumn = NULL, idType = NULL, pvalColumn = NULL, logPvalColumn = NULL, log2FCColumn = NULL, signalColumn = NULL )
getMetDEMeta( met.de.raw, met.db, idColumn = NULL, idType = NULL, pvalColumn = NULL, logPvalColumn = NULL, log2FCColumn = NULL, signalColumn = NULL )
met.de.raw |
A table with differential expression results, an object convertable to data.frame. |
met.db |
Metabolite database |
idColumn |
Specifies column name with metabolite identifiers. |
idType |
Specifies type of metabolite IDs (one of the supported by annotation). |
pvalColumn |
Specifies column with p-values. |
logPvalColumn |
Specifies column with log p-values, if there is no such column one will be generated automatically. |
log2FCColumn |
Specifies column with log2-fold changes. |
signalColumn |
Specifies column with identifier of the measured entity Could be NULL (automatic, set from based on pval and log2FC columns), character (column name), or function (evaluated in a scope of original data frame) |
object with prepared columns for the analysis for metabolite data
data("met.kegg.dbEx") data("met.de.rawEx") de.meta <- getMetDEMeta(met.de.rawEx, met.db = met.kegg.dbEx)
data("met.kegg.dbEx") data("met.de.rawEx") de.meta <- getMetDEMeta(met.de.rawEx, met.db = met.kegg.dbEx)
See file https://github.com/ctlab/gatom/blob/master/inst/scripts/example.R for details.
igraph object
See file https://github.com/ctlab/gatom/blob/master/inst/scripts/example.R for details.
igraph object
Creates metabolic graph based on specified data
makeMetabolicGraph( network, topology = c("atoms", "metabolites"), org.gatom.anno, gene.de, gene.de.meta = getGeneDEMeta(gene.de, org.gatom.anno), gene.keep.top = 12000, met.db, met.de, met.de.meta = getMetDEMeta(met.de, met.db), met.to.filter = fread(system.file("extdata", "mets2mask.lst", package = "gatom"))$ID, gene2reaction.extra = NULL, keepReactionsWithoutEnzymes = FALSE, largest.component = TRUE )
makeMetabolicGraph( network, topology = c("atoms", "metabolites"), org.gatom.anno, gene.de, gene.de.meta = getGeneDEMeta(gene.de, org.gatom.anno), gene.keep.top = 12000, met.db, met.de, met.de.meta = getMetDEMeta(met.de, met.db), met.to.filter = fread(system.file("extdata", "mets2mask.lst", package = "gatom"))$ID, gene2reaction.extra = NULL, keepReactionsWithoutEnzymes = FALSE, largest.component = TRUE )
network |
Network object |
topology |
Way to determine network vertices |
org.gatom.anno |
Organism annotation object |
gene.de |
Table with the differential gene expression, set to NULL if absent |
gene.de.meta |
Annotation of 'gene.de' table |
gene.keep.top |
Only the 'gene.keep.top' of the most expressed genes will be kept for the network |
met.db |
Metabolite database |
met.de |
Table with the differential expression for metabolites, set to NULL if absent |
met.de.meta |
Annotation of 'met.de' table |
met.to.filter |
List of metabolites to filter from the network |
gene2reaction.extra |
Additional gene to reaction mappings. Should be a data.table with 'gene' and 'reaction' columns |
keepReactionsWithoutEnzymes |
If TRUE, keep reactions that have no annotated enzymes, thus expanding the network but including some reactions which are not possible in the considered species. |
largest.component |
If TRUE, only the largest connected component is returned |
igraph object created from input data
data("gene.de.rawEx") data("met.de.rawEx") data("met.kegg.dbEx") data("networkEx") data("org.Mm.eg.gatom.annoEx") g <- makeMetabolicGraph(network = networkEx, topology = "atoms", org.gatom.anno = org.Mm.eg.gatom.annoEx, gene.de = gene.de.rawEx, met.db = met.kegg.dbEx, met.de = met.de.rawEx)
data("gene.de.rawEx") data("met.de.rawEx") data("met.kegg.dbEx") data("networkEx") data("org.Mm.eg.gatom.annoEx") g <- makeMetabolicGraph(network = networkEx, topology = "atoms", org.gatom.anno = org.Mm.eg.gatom.annoEx, gene.de = gene.de.rawEx, met.db = met.kegg.dbEx, met.de = met.de.rawEx)
Create an organism annotation object for network analysis
makeOrgGatomAnnotation( org.db, idColumns = c(Entrez = "ENTREZID", RefSeq = "REFSEQ", Ensembl = "ENSEMBL", Symbol = "SYMBOL"), nameColumn = "SYMBOL", enzymeColumn = "ENZYME", appendEnzymesFromKegg = TRUE, appendOrthologiesFromKegg = TRUE, filterNonSpecificEnzymes = TRUE, keggOrgCode = NULL )
makeOrgGatomAnnotation( org.db, idColumns = c(Entrez = "ENTREZID", RefSeq = "REFSEQ", Ensembl = "ENSEMBL", Symbol = "SYMBOL"), nameColumn = "SYMBOL", enzymeColumn = "ENZYME", appendEnzymesFromKegg = TRUE, appendOrthologiesFromKegg = TRUE, filterNonSpecificEnzymes = TRUE, keggOrgCode = NULL )
org.db |
Bioconductor org.db object, e.g. org.Mm.eg.db |
idColumns |
vector of column names from 'org.db' object to creat ID mappings. First ID will be used as a base identifier, should be compatible with KEGG and Reactome databases. |
nameColumn |
column with a human readable gene symbol. Default to "SYMBOL". |
enzymeColumn |
column with an Enzyme Commission ID. Default to "ENZYME". |
appendEnzymesFromKegg |
if TRUE, KEGG databases will be sued to extend gene to enzyme mappings obtained from org.db package. |
appendOrthologiesFromKegg |
if TRUE, KEGG database will be sued to extend gene to orthology mappings obtained from org.db package |
filterNonSpecificEnzymes |
if TRUE, will filter out non-specific enzymes from gene to enzyme mappings obtained from org.db package |
keggOrgCode |
KEGG organism code, e.g. "mmu". If set to NULL, the code is determined automatically. |
organism annotation object that will be used for network analysis
library(org.Mm.eg.db) org.Mm.eg.gatom.anno <- makeOrgGatomAnnotation(org.db = org.Mm.eg.db)
library(org.Mm.eg.db) org.Mm.eg.gatom.anno <- makeOrgGatomAnnotation(org.db = org.Mm.eg.db)
See file https://github.com/ctlab/gatom/blob/master/inst/scripts/example.R for details.
tibble/data.frame object
See file https://github.com/ctlab/gatom/blob/master/inst/scripts/example.R for details.
list object
See file https://github.com/ctlab/gatom/blob/master/inst/scripts/example.R for details.
igraph object
See file https://github.com/ctlab/gatom/blob/master/inst/scripts/example.R for details.
list object
See file https://github.com/ctlab/gatom/blob/master/inst/scripts/example.R for details.
list object
Makes data.table with differential expression results containing all columns required for gatom and in the expected format based on metadata object
prepareDE(de.raw, de.meta)
prepareDE(de.raw, de.meta)
de.raw |
Table with differential expression results, an object convertable to data.frame |
de.meta |
Object with differential expression table metadata acquired with getGeneDEMeta or getMetDEMeta functions |
data.table object with converted differential expression table
data("org.Mm.eg.gatom.annoEx") data("gene.de.rawEx") de.meta <- getGeneDEMeta(gene.de.rawEx, org.gatom.anno = org.Mm.eg.gatom.annoEx) de <- prepareDE(gene.de.rawEx, de.meta)
data("org.Mm.eg.gatom.annoEx") data("gene.de.rawEx") de.meta <- getGeneDEMeta(gene.de.rawEx, org.gatom.anno = org.Mm.eg.gatom.annoEx) de <- prepareDE(gene.de.rawEx, de.meta)
Save module to a graphviz dot file
saveModuleToDot( module, file, name = NULL, extra.node.attrs = NULL, extra.edge.attrs = NULL )
saveModuleToDot( module, file, name = NULL, extra.node.attrs = NULL, extra.edge.attrs = NULL )
module |
Module to save |
file |
File to save to |
name |
Name of the module |
extra.node.attrs |
Table with additional node attributes to be written to the dot file as is |
extra.edge.attrs |
Table with additional edge attributes to be written to the dot file as is |
Returns NULL
data(mEx) saveModuleToDot(module = mEx, file = "module.dot")
data(mEx) saveModuleToDot(module = mEx, file = "module.dot")
Save module to a html widget
saveModuleToHtml( module, file, name = "", sizingPolicy = htmlwidgets::sizingPolicy(defaultWidth = "100%", defaultHeight = "90vh", padding = 10), ... )
saveModuleToHtml( module, file, name = "", sizingPolicy = htmlwidgets::sizingPolicy(defaultWidth = "100%", defaultHeight = "90vh", padding = 10), ... )
module |
Module to save |
file |
File to save to |
name |
Name of the module |
sizingPolicy |
A widget sizing policy |
... |
Other parameters |
Returns NULL
data(mEx) saveModuleToHtml(module = mEx, file = "module.html")
data(mEx) saveModuleToHtml(module = mEx, file = "module.html")
Save module to a nice pdf file
saveModuleToPdf(module, file, name = NULL, n_iter = 100, force = 1e-05)
saveModuleToPdf(module, file, name = NULL, n_iter = 100, force = 1e-05)
module |
Module to save |
file |
File to save to |
name |
Name of the module |
n_iter |
Number of repel algorithm iterations |
force |
Value of repel force |
Returns NULL
data(mEx) saveModuleToPdf(module = mEx, file = "module.pdf")
data(mEx) saveModuleToPdf(module = mEx, file = "module.pdf")
Save module to an XGMML file
saveModuleToXgmml(module, file, name = NULL)
saveModuleToXgmml(module, file, name = NULL)
module |
Module to save |
file |
File to save to |
name |
Name of the module |
Returns NULL
data(mEx) saveModuleToXgmml(module = mEx, file = "module.xgmml")
data(mEx) saveModuleToXgmml(module = mEx, file = "module.xgmml")
Score metabolic graph
scoreGraph( g, k.gene, k.met, vertex.threshold.min = 0.1, edge.threshold.min = 0.1, met.score.coef = 1, show.warnings = TRUE, raw = FALSE )
scoreGraph( g, k.gene, k.met, vertex.threshold.min = 0.1, edge.threshold.min = 0.1, met.score.coef = 1, show.warnings = TRUE, raw = FALSE )
g |
Metabolic graph obtained with makeMetabolic graph function |
k.gene |
Number of gene signals to be scored positively, the higher is the number, the larger will be the resulting module. If set to NULL, genes will not be used for scoring. |
k.met |
Number of metabolite signals to be scored positively, the higher is the number, the larger will be the resulting module. If set to NULL, metabolites will not be used for scoring. |
vertex.threshold.min |
The worst acceptable estimated FDR for vertices. If necessary number of positive metabolite signals will be decreased from 'k.met' to reach this threshold. Default value is 0.1. |
edge.threshold.min |
The worst acceptable estimated FDR for vertices. If necessary number of positive metabolite signals will be decreased from 'k.gene' to reach this threshold. Default value is 0.1. |
met.score.coef |
Coefficient on which all vertex weights are multiplied. Can be used to balance vertex and edge weights. Default values is 1. |
show.warnings |
whether to show warnings |
raw |
whether to return raw scored graph, not a SGMWCS instance. Default to FALSE. |
SGMWCS instance or scored igraph object
data("gEx") gs <- scoreGraph(g = gEx, k.gene = 25, k.met = 25)
data("gEx") gs <- scoreGraph(g = gEx, k.gene = 25, k.met = 25)