Title: | NetPathMiner for Biological Network Construction, Path Mining and Visualization |
---|---|
Description: | NetPathMiner is a general framework for network path mining using genome-scale networks. It constructs networks from KGML, SBML and BioPAX files, providing three network representations, metabolic, reaction and gene representations. NetPathMiner finds active paths and applies machine learning methods to summarize found paths for easy interpretation. It also provides static and interactive visualizations of networks and paths to aid manual investigation. |
Authors: | Ahmed Mohamed [aut, cre] , Tim Hancock [aut], Tim Hancock [aut] |
Maintainer: | Ahmed Mohamed <[email protected]> |
License: | GPL (>= 2) |
Version: | 1.43.0 |
Built: | 2024-10-30 08:27:04 UTC |
Source: | https://github.com/bioc/NetPathMiner |
NetPathMiner implements a flexible module-based process flow for network path mining and visualization, which can be fully inte-grated with user-customized functions. NetPathMiner supports construction of various types of genome scale networks from KGML, SBML and BioPAX formats, enabling its utility to most common pathway databases. NetPathMiner also provides different visualization techniques to facilitate the analysis of even thousands of output paths.
Ahmed Mohamed [email protected]
This function computes edge weights based on a gene expression profile.
assignEdgeWeights( microarray, graph, use.attr, y, weight.method = "cor", complex.method = "max", missing.method = "median", same.gene.penalty = "median", bootstrap = 100, verbose = TRUE )
assignEdgeWeights( microarray, graph, use.attr, y, weight.method = "cor", complex.method = "max", missing.method = "median", same.gene.penalty = "median", bootstrap = 100, verbose = TRUE )
microarray |
Microarray should be a Dataframe or a matrix, with genes as rownames, and samples as columns. |
graph |
An annotated igraph object. |
use.attr |
An attribute name to map |
y |
Sample labels, given as a factor or a character vector. This must be the same size as the columns of |
weight.method |
A function, or a string indicating the name of the function to be used to compute the edge weights.
The function is provided with 2 numerical verctors (2 rows from |
complex.method |
A function, or a string indicating the name of the function to be used in weighting edges connecting complexes.
If a vertex has >1 attribute value, all possible pairwise weights are first computed, and given to |
missing.method |
A function, or a string indicating the name of the function to be used in weighting edges when one of the vertices
lack expression data. The function is passed all edge weights on the graph. Default is |
same.gene.penalty |
A numerical value to be assigned when 2 adjacent vertices have the same attribute value, since correlation and
similarity measure will give perfect scores. Alternatively, |
bootstrap |
An integer |
verbose |
Print the progress of the function. |
The input graph with edge.weight
as an edge attribute. The attribute can be a list of weights if y
labels
were provided.
Ahmed Mohamed
## Convert a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) # Using Spearman correlation, assigning missing edges to -1 ## Not run: assignEdgeWeights(microarray, graph, use.attr="miriam.affy.probeset", y=factor(colnames(microarray)), weight.method = function(x1,x2) cor(x1,x2, method="spearman"), missing.method = -1) ## End(Not run)
## Convert a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) # Using Spearman correlation, assigning missing edges to -1 ## Not run: assignEdgeWeights(microarray, graph, use.attr="miriam.affy.probeset", y=factor(colnames(microarray)), weight.method = function(x1,x2) cor(x1,x2, method="spearman"), missing.method = -1) ## End(Not run)
This function takes BioPAX objects (level 2 or 3) as input, and returns either a metabolic or a signaling network as output.
biopax2igraph( biopax, parse.as = c("metabolic", "signaling"), expand.complexes = FALSE, inc.sm.molecules = FALSE, verbose = TRUE )
biopax2igraph( biopax, parse.as = c("metabolic", "signaling"), expand.complexes = FALSE, inc.sm.molecules = FALSE, verbose = TRUE )
biopax |
BioPAX object generated by |
parse.as |
Whether to process file into a metabolic or a signaling network. |
expand.complexes |
Split protein complexes into individual gene nodes. Ignored if
|
inc.sm.molecules |
Include small molecules that are participating in signaling events. Ignored if
|
verbose |
Whether to display the progress of the function. |
This function requires rBiopaxParser
installed.
Users can specify whether files are processes as metabolic or signaling networks.
Metabolic networks are given as bipartite graphs, where metabolites and reactions represent
vertex types. Reactions are constructed from Conversion
classes, connecting them
to their corresponding Left
s and Right
s. Each reaction vertex has genes
attribute,
listing all Catalysis
relationships of this reaction. As a general rule, reactions inherit all annotation
attributes of its catalyzig genes.
Signaling network have genes as vertices and edges represent interactions, such as activiation / inhibition.
Genes participating in successive reactions are also connected. Signaling interactions are constructed from
Control
classes, where edges are drawn from controller
to controlled
.
All annotation attributes are exracted from XRefs
associated with the vertices, and are stored according to
MIRIAM guidelines (miraim.db
, where db is the database name).
An igraph object, representing a metbolic or a signaling network.
Ahmed Mohamed
Other Database extraction methods:
KGML2igraph()
,
SBML2igraph()
if(requireNamespace("rBiopaxParser")){ data(ex_biopax) # Process biopax as a metabolic network g <- biopax2igraph(ex_biopax) plotNetwork(g) # Process SBML file as a signaling network g <- biopax2igraph(ex_biopax, parse.as="signaling", expand.complexes=TRUE) }
if(requireNamespace("rBiopaxParser")){ data(ex_biopax) # Process biopax as a metabolic network g <- biopax2igraph(ex_biopax) plotNetwork(g) # Process SBML file as a signaling network g <- biopax2igraph(ex_biopax, parse.as="signaling", expand.complexes=TRUE) }
This function returns a list of colors for vertices, assigned similar colors if they share a common attribute (ex: in the same pathway, etc).
colorVertexByAttr(graph, attr.name, col.palette = palette())
colorVertexByAttr(graph, attr.name, col.palette = palette())
graph |
An annotated igraph object. |
attr.name |
The attribute name (ex: "pathway") by which vertices will be colored. Complex attributes, where a vertex belongs to more than one group, are supported. |
col.palette |
A color palette, or a palette generating function (ex: col.palette=rainbow ). |
A list of colors (in HEX format) for vertices.
Ahmed Mohamed
Other Plotting methods:
layoutVertexByAttr()
,
plotAllNetworks()
,
plotClassifierROC()
,
plotClusterMatrix()
,
plotCytoscapeGML()
,
plotNetwork()
,
plotPathClassifier()
,
plotPaths()
data("ex_kgml_sig") v.colors <- colorVertexByAttr(ex_kgml_sig, "pathway") plotNetwork(ex_kgml_sig, vertex.color=v.colors)
data("ex_kgml_sig") v.colors <- colorVertexByAttr(ex_kgml_sig, "pathway") plotNetwork(ex_kgml_sig, vertex.color=v.colors)
A dataset containing Porphyrin metabolism pathway in Biopax Level 3 and parsed with
readBiopax
.
data(ex_biopax) ex_biopax
data(ex_biopax) ex_biopax
An example igraph object representing Ras and chemokine signaling pathways in human extracted from KGML files.
data(ex_kgml_sig) plotNetwork(ex_kgml_sig, vertex.color="pathway")
data(ex_kgml_sig) plotNetwork(ex_kgml_sig, vertex.color="pathway")
An microarray data example. This is part of the ALL dataset, for demonstration purposes.
data(ex_microarray)
data(ex_microarray)
An example igraph object representing bipartite metabolic network of Carbohydrate metabolism extracted from SBML file from Reactome database.
data(ex_sbml) plotNetwork(ex_sbml, vertex.color="compartment.name")
data(ex_sbml) plotNetwork(ex_sbml, vertex.color="compartment.name")
These are general functions to expand vertices by their attributes, i.e. create a separate vertex for each attribute value.
expandComplexes( graph, v.attr, keep.parent.attr = "^pathway", expansion.method = c("normal", "duplicate"), missing.method = c("keep", "remove", "reconnect") ) makeGeneNetwork( graph, v.attr = "genes", keep.parent.attr = "^pathway", expansion.method = "duplicate", missing.method = "remove" )
expandComplexes( graph, v.attr, keep.parent.attr = "^pathway", expansion.method = c("normal", "duplicate"), missing.method = c("keep", "remove", "reconnect") ) makeGeneNetwork( graph, v.attr = "genes", keep.parent.attr = "^pathway", expansion.method = "duplicate", missing.method = "remove" )
graph |
An annotated igraph object. |
v.attr |
Name of the attribute which vertices are expanded to. |
keep.parent.attr |
A (List of) |
expansion.method |
If |
missing.method |
How to deal with vertices with no attribute values. |
These functions can be very useful when merging networks constructed from different databases. For example, to match a network created from Reactome to a KEGG network, you can expand metabolite vertices by "miriam.kegg.compound" attribute.
A new graph with vertices expanded.
makeGeneNetwork
returns a graph, where nodes are genes, and edges represent
participation in succesive reactions.
Ahmed Mohamed
Other Network processing methods:
makeMetaboliteNetwork()
,
makeReactionNetwork()
,
reindexNetwork()
,
rmSmallCompounds()
,
simplifyReactionNetwork()
,
vertexDeleteReconnect()
## Make a gene network from a reaction network. data(ex_sbml) # A bipartite metbaolic network. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ggraph <- makeGeneNetwork(rgraph) ## Expand vertices into their contituent genes. data(ex_kgml_sig) # Ras and chemokine signaling pathways in human ggraph <- expandComplexes(ex_kgml_sig, v.attr = "miriam.ncbigene", keep.parent.attr= c("^pathway", "^compartment")) ## Create a separate vertex for each compartment. This is useful in duplicating ## metabolite vertices in a network. ## Not run: graph <- expandComplexes(graph, v.attr = "compartment", keep.parent.attr = "all", expansion.method = "duplicate", missing.method = "keep") ## End(Not run)
## Make a gene network from a reaction network. data(ex_sbml) # A bipartite metbaolic network. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ggraph <- makeGeneNetwork(rgraph) ## Expand vertices into their contituent genes. data(ex_kgml_sig) # Ras and chemokine signaling pathways in human ggraph <- expandComplexes(ex_kgml_sig, v.attr = "miriam.ncbigene", keep.parent.attr= c("^pathway", "^compartment")) ## Create a separate vertex for each compartment. This is useful in duplicating ## metabolite vertices in a network. ## Not run: graph <- expandComplexes(graph, v.attr = "compartment", keep.parent.attr = "all", expansion.method = "duplicate", missing.method = "keep") ## End(Not run)
Creates a subnetwork from a ranked path list generated by pathRanker
.
extractPathNetwork(paths, graph)
extractPathNetwork(paths, graph)
paths |
The paths extracted by |
graph |
A annotated igraph object. |
A subnetwork from all paths provided. If paths are computed for several labels (sample categories), a subnetwork is returned for each label.
Ahmed Mohamed
Other Path ranking methods:
getPathsAsEIDs()
,
pathRanker()
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=6) ## Get the subnetwork of paths in reaction graph. reaction.sub <- getPathsAsEIDs(ranked.p, rgraph) ## Get the subnetwork of paths in the original metabolic graph. metabolic.sub <- getPathsAsEIDs(ranked.p, ex_sbml)
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=6) ## Get the subnetwork of paths in reaction graph. reaction.sub <- getPathsAsEIDs(ranked.p, rgraph) ## Get the subnetwork of paths in the original metabolic graph. metabolic.sub <- getPathsAsEIDs(ranked.p, ex_sbml)
These functions report the annotation status of the vertices of a given network, modify or remove certain annotations.
getAttrStatus(graph, pattern = "^miriam.") getAttrNames(graph, pattern = "") getAttribute(graph, attr.name) setAttribute(graph, attr.name, attr.value) rmAttribute(graph, attr.name)
getAttrStatus(graph, pattern = "^miriam.") getAttrNames(graph, pattern = "") getAttribute(graph, attr.name) setAttribute(graph, attr.name, attr.value) rmAttribute(graph, attr.name)
graph |
An annotated igraph object. |
pattern |
A |
attr.name |
The attribute name |
attr.value |
A list of attribute values. This must be the same size as the number of vertices. |
NetPathMiner stores all its vertex annotation attributes in a list, and stores them collectively as
a single attr
. This is not to interfer with graph_attr_names
from igraph
package.
All functions here target NetPathMiner annotations only.
For getAttrStatus
, a dataframe summarizing the number of vertices with no (missing
), one (single
)
or more than one (complex
) attribute value. The coverage
For getAttrNames
, a character vector of attribute names matching the pattern.
For getAttribute
, a list of vertex annotation values for the query attribute.
For setAttribute
, a graph with the new attribute set.
For rmAttrNames
, a new igraph object with the attibute removed.
Ahmed Mohamed
Other Attribute handling methods:
stdAttrNames()
data(ex_kgml_sig) # Ras and chemokine signaling pathways in human # Get status of attribute "pathway" only getAttrStatus(ex_kgml_sig, "^pathway$") # Get status of all attributes starting with "pathway" and "miriam" keywords getAttrStatus(ex_kgml_sig, "(^miriam)|(^pathway)") # Get all attribute names containing "miriam" getAttrNames(ex_kgml_sig, "miriam") # Get all attribute names containing "miriam" getAttribute(ex_kgml_sig, "miriam.ncbigene") # Remove an attribute from graph graph <- rmAttribute(ex_kgml_sig, "miriam.ncbigene")
data(ex_kgml_sig) # Ras and chemokine signaling pathways in human # Get status of attribute "pathway" only getAttrStatus(ex_kgml_sig, "^pathway$") # Get status of all attributes starting with "pathway" and "miriam" keywords getAttrStatus(ex_kgml_sig, "(^miriam)|(^pathway)") # Get all attribute names containing "miriam" getAttrNames(ex_kgml_sig, "miriam") # Get all attribute names containing "miriam" getAttribute(ex_kgml_sig, "miriam.ncbigene") # Remove an attribute from graph graph <- rmAttribute(ex_kgml_sig, "miriam.ncbigene")
This function generates geneset networks based on a given netowrk, by grouping vertices sharing common attributes (in the same pathway or compartment).
getGeneSetNetworks( graph, use.attr = "pathway", format = c("list", "pathway-class") )
getGeneSetNetworks( graph, use.attr = "pathway", format = c("list", "pathway-class") )
graph |
An annotated igraph object.. |
use.attr |
The attribute by which vertices are grouped (tepically pathway, or GO) |
format |
The output format. If "list" is specified, a list of subgraphs are returned (default). If "pathway-class" is specified, a list of pathway-class objects are returned. Pathway-class is used by graphite package to run several methods of topology-based enrichment analyses. |
A list of geneset networks as igraph or Pathway-class objects.
Ahmed Mohamed
data(ex_kgml_sig) # Ras and chemokine signaling pathways in human genesetnets <- getGeneSetNetworks(ex_kgml_sig, use.attr="pathway") # Integration with graphite package ## Not run: if(requireNamespace("graphite") & requireNamespace("clipper") & requireNamespace("ALL")){ genesetnets <- getGeneSetNetworks(ex_kgml_sig, use.attr="pathway", format="pathway-class") path <- convertIdentifiers(genesetnets$`Chemokine signaling pathway`, "entrez") genes <- nodes(path) data(ALL) all <- as.matrix(exprs(ALL[1:length(genes),1:20])) classes <- c(rep(1,10), rep(2,10)) rownames(all) <- genes runClipper(path, all, classes, "mean", pathThr=0.1) } ## End(Not run)
data(ex_kgml_sig) # Ras and chemokine signaling pathways in human genesetnets <- getGeneSetNetworks(ex_kgml_sig, use.attr="pathway") # Integration with graphite package ## Not run: if(requireNamespace("graphite") & requireNamespace("clipper") & requireNamespace("ALL")){ genesetnets <- getGeneSetNetworks(ex_kgml_sig, use.attr="pathway", format="pathway-class") path <- convertIdentifiers(genesetnets$`Chemokine signaling pathway`, "entrez") genes <- nodes(path) data(ALL) all <- as.matrix(exprs(ALL[1:length(genes),1:20])) classes <- c(rep(1,10), rep(2,10)) rownames(all) <- genes runClipper(path, all, classes, "mean", pathThr=0.1) } ## End(Not run)
This function generates genesets based on a given netowrk, by grouping vertices sharing
common attributes (in the same pathway or compartment). Genes associated with each vertex
can be specified through gene.attr
argument.
getGeneSets(graph, use.attr = "pathway", gene.attr = "genes", gmt.file)
getGeneSets(graph, use.attr = "pathway", gene.attr = "genes", gmt.file)
graph |
An annotated igraph object.. |
use.attr |
The attribute by which vertices are grouped (tepically pathway, or GO) |
gene.attr |
The attribute listing genes annotated with each vertex (ex: miriam.ncbigene, miriam.uniprot, ...) |
gmt.file |
Optinal. If provided, Results are exported to a GMT file. GMT files are readily used by most gene set analysis packages. |
A list of genesets or written to gmt file if provided.
Ahmed Mohamed
data(ex_kgml_sig) # Ras and chemokine signaling pathways in human genesets <- getGeneSets(ex_kgml_sig, use.attr="pathway", gene.attr="miriam.ncbigene") # Write the genesets in a GMT file, and read it using GSEABase package. getGeneSets(ex_kgml_sig, use.attr="pathway", gene.attr="miriam.ncbigene", gmt.file="kgml.gmt") ## Not run: if(requireNamespace("GSEABase")) toGmt("kgml.gmt") ## End(Not run) # Create genesets using compartment information data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. genesets <- getGeneSets(ex_sbml, use.attr="compartment.name", gene.attr="miriam.uniprot")
data(ex_kgml_sig) # Ras and chemokine signaling pathways in human genesets <- getGeneSets(ex_kgml_sig, use.attr="pathway", gene.attr="miriam.ncbigene") # Write the genesets in a GMT file, and read it using GSEABase package. getGeneSets(ex_kgml_sig, use.attr="pathway", gene.attr="miriam.ncbigene", gmt.file="kgml.gmt") ## Not run: if(requireNamespace("GSEABase")) toGmt("kgml.gmt") ## End(Not run) # Create genesets using compartment information data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. genesets <- getGeneSets(ex_sbml, use.attr="compartment.name", gene.attr="miriam.uniprot")
Convert a ranked path list to Edge ids of a graph, where paths can come from a different representation (for example matching path from a reaction network to edges on a metabolic network).
getPathsAsEIDs(paths, graph)
getPathsAsEIDs(paths, graph)
paths |
The paths extracted by |
graph |
A annotated igraph object. |
A list of edge ids on the provided graph.
Ahmed Mohamed
Other Path ranking methods:
extractPathNetwork()
,
pathRanker()
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=6) ## Get the edge ids along paths in the reaction graph. path.eids <- getPathsAsEIDs(ranked.p, rgraph) ## Get the edge ids along paths in the original metabolic graph. path.eids <- getPathsAsEIDs(ranked.p, ex_sbml)
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=6) ## Get the edge ids along paths in the reaction graph. path.eids <- getPathsAsEIDs(ranked.p, rgraph) ## Get the edge ids along paths in the original metabolic graph. path.eids <- getPathsAsEIDs(ranked.p, ex_sbml)
This function takes KGML files as input, and returns either a metabolic or a signaling network as output.
KGML2igraph( filename, parse.as = c("metabolic", "signaling"), expand.complexes = FALSE, verbose = TRUE )
KGML2igraph( filename, parse.as = c("metabolic", "signaling"), expand.complexes = FALSE, verbose = TRUE )
filename |
A character vector containing the KGML files to be processed. If a directory path is provided, all *.xml files in it and its subdirectories are included. |
parse.as |
Whether to process file into a metabolic or a signaling network. |
expand.complexes |
Split protein complexes into individual gene nodes. This argument is
ignored if |
verbose |
Whether to display the progress of the function. |
Users can specify whether files are processes as metabolic or signaling networks.
Metabolic networks are given as bipartite graphs, where metabolites and reactions represent
vertex types. This is constructed from <reaction> xml node in KGML file, connecting them
to their corresponding substrates and products. Each reaction vertex has genes
attribute,
listing all genes associated with the reaction. As a general rule, reactions inherit all annotation
attributes of its catalyzig genes.
Signaling network have genes as vertices and edges represent interactions, such as activiation / inhibition. Genes participating in successive reactions are also connected. Signaling parsing method processes <ECrel>, <PPrel> and <PCrel> interactions from KGML files.
To generate a genome scale network, simply provide a list of files to be parsed, or put all
file in a directory, as pass the directory path as filename
An igraph object, representing a metbolic or a signaling network.
Ahmed Mohamed
Other Database extraction methods:
SBML2igraph()
,
biopax2igraph()
if(is.loaded("readkgmlfile")){ # This is false if libxml2 wasn't available at installation. filename <- system.file("extdata", "hsa00860.xml", package="NetPathMiner") # Process KGML file as a metabolic network g <- KGML2igraph(filename) plotNetwork(g) # Process KGML file as a signaling network g <- KGML2igraph(filename, parse.as="signaling", expand.complexes=TRUE) plotNetwork(g) }
if(is.loaded("readkgmlfile")){ # This is false if libxml2 wasn't available at installation. filename <- system.file("extdata", "hsa00860.xml", package="NetPathMiner") # Process KGML file as a metabolic network g <- KGML2igraph(filename) plotNetwork(g) # Process KGML file as a signaling network g <- KGML2igraph(filename, parse.as="signaling", expand.complexes=TRUE) plotNetwork(g) }
This function generates a layout for igraph objects, keeping vertices with the same attribute (ex: in the same pathway, etc) close to each other.
layoutVertexByAttr( graph, attr.name, cluster.strength = 1, layout = layout.auto )
layoutVertexByAttr( graph, attr.name, cluster.strength = 1, layout = layout.auto )
graph |
An annotated igraph object. |
attr.name |
The attribute name by which vertices are laid out. |
cluster.strength |
A number indicating tie strengths between vertices with the same attribute. The larger it is, the closer the vertices will be. |
layout |
A layout function, ideally a force-directed layout fuction, such as
|
A two-column matrix indicating the x and y postions of vertices.
Ahmed Mohamed
Other Plotting methods:
colorVertexByAttr()
,
plotAllNetworks()
,
plotClassifierROC()
,
plotClusterMatrix()
,
plotCytoscapeGML()
,
plotNetwork()
,
plotPathClassifier()
,
plotPaths()
data("ex_kgml_sig") v.layout <- layoutVertexByAttr(ex_kgml_sig, "pathway") plotNetwork(ex_kgml_sig, vertex.color="pathway", layout=v.layout) v.layout <- layoutVertexByAttr(ex_kgml_sig, "pathway", cluster.strength=5) plotNetwork(ex_kgml_sig, vertex.color="pathway", layout=v.layout)
data("ex_kgml_sig") v.layout <- layoutVertexByAttr(ex_kgml_sig, "pathway") plotNetwork(ex_kgml_sig, vertex.color="pathway", layout=v.layout) v.layout <- layoutVertexByAttr(ex_kgml_sig, "pathway", cluster.strength=5) plotNetwork(ex_kgml_sig, vertex.color="pathway", layout=v.layout)
This function removes reaction nodes keeping them as edge attributes. The resulting network contains metabolite nodes only, where edges indicate that reaction conversions.
makeMetaboliteNetwork(graph)
makeMetaboliteNetwork(graph)
graph |
A metabolic network. |
A reaction network.
Ahmed Mohamed
Other Network processing methods:
expandComplexes()
,
makeReactionNetwork()
,
reindexNetwork()
,
rmSmallCompounds()
,
simplifyReactionNetwork()
,
vertexDeleteReconnect()
## Conver a metabolic network to a metbolite network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. mgraph <- makeMetaboliteNetwork(ex_sbml)
## Conver a metabolic network to a metbolite network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. mgraph <- makeMetaboliteNetwork(ex_sbml)
This function removes metabolite nodes keeping them as edge attributes. The resulting network contains reaction nodes only, where edges indicate that a metabolite produced by one reaction is consumed by the other.
makeReactionNetwork(graph, simplify = FALSE)
makeReactionNetwork(graph, simplify = FALSE)
graph |
A metabolic network. |
simplify |
An option to remove translocation and spontaneous reactions that require no catalyzing genes. Translocation reactions are detected from reaction name (SBML, BioPAX), or by having identical substrates and products. |
A reaction network.
Ahmed Mohamed
Other Network processing methods:
expandComplexes()
,
makeMetaboliteNetwork()
,
reindexNetwork()
,
rmSmallCompounds()
,
simplifyReactionNetwork()
,
vertexDeleteReconnect()
## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE)
## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE)
This function gets a NetPathMiner default value for a variable.
NPMdefaults(value)
NPMdefaults(value)
value |
a character string indicating the variable name. |
NetPathMiner defines the following defaults:
small.comp.ls Dataframe of ubiquitous metabolites. Used by rmSmallCompounds
.
bridge Dataframe of attributes supported by Brigde Database. Used by fetchAttribute
.
bridge.organisms A list of bridge supported organisms. Used by fetchAttribute
.
bridge.web The base URL for Brigde Database webservices. Used by fetchAttribute
.
The defuult value for the given variable.
Ahmed Mohamed
# Get the default list of small compounds (uniquitous metabolites). NPMdefaults("small.comp.ls")
# Get the default list of small compounds (uniquitous metabolites). NPMdefaults("small.comp.ls")
HME3M Markov pathway classifier.
pathClassifier( paths, target.class, M, alpha = 1, lambda = 2, hme3miter = 100, plriter = 1, init = "random" )
pathClassifier( paths, target.class, M, alpha = 1, lambda = 2, hme3miter = 100, plriter = 1, init = "random" )
paths |
The training paths computed by |
target.class |
he label of the targe class to be classified. This label must be present
as a label within the |
M |
Number of components within the paths to be extracted. |
alpha |
The PLR learning rate. (between 0 and 1). |
lambda |
The PLR regularization parameter. (between 0 and 2) |
hme3miter |
Maximum number of HME3M iterations. It will stop when likelihood change is < 0.001. |
plriter |
Maximum number of PLR iteractions. It will stop when likelihood change is < 0.001. |
init |
Specify whether to initialize the HME3M responsibilities with the 3M model - random is recommended. |
Take care with selection of lambda and alpha - make sure you check that the likelihood is always increasing.
A list with the following elements. A list with the following values
h |
A dataframe with the EM responsibilities. |
theta |
A dataframe with the Markov parameters for each component. |
beta |
A dataframe with the PLR coefficients for each component. |
proportions |
The probability of each HME3M component. |
posterior.probs |
The HME3M posterior probability. |
likelihood |
The likelihood convergence history. |
plrplr |
The posterior predictions from each components PLR model. |
path.probabilities |
The 3M probabilities for each path belonging to each component. |
params |
The parameters used to build the model. |
y |
The binary response variable used by HME3M. A 1 indicates the location of the target.class labels in |
perf |
The training set ROC curve AUC. |
label |
The HME3M predicted label for each path. |
component |
The HME3M component assignment for each path. |
Timothy Hancock and Ichigaku Takigawa
Hancock, Timothy, and Mamitsuka, Hiroshi: A Markov Classification Model for Metabolic Pathways, Workshop on Algorithms in Bioinformatics (WABI) , 2009
Hancock, Timothy, and Mamitsuka, Hiroshi: A Markov Classification Model for Metabolic Pathways, Algorithms for Molecular Biology 2010
Other Path clustering & classification methods:
pathCluster()
,
pathsToBinary()
,
plotClassifierROC()
,
plotClusterMatrix()
,
plotPathClassifier()
,
plotPathCluster()
,
predictPathClassifier()
,
predictPathCluster()
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=6) ## Convert paths to binary matrix. ybinpaths <- pathsToBinary(ranked.p) p.class <- pathClassifier(ybinpaths, target.class = "BCR/ABL", M = 3) ## Contingency table of classification performance table(ybinpaths$y,p.class$label) ## Plotting the classifier results. plotClassifierROC(p.class) plotClusters(ybinpaths, p.class)
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=6) ## Convert paths to binary matrix. ybinpaths <- pathsToBinary(ranked.p) p.class <- pathClassifier(ybinpaths, target.class = "BCR/ABL", M = 3) ## Contingency table of classification performance table(ybinpaths$y,p.class$label) ## Plotting the classifier results. plotClassifierROC(p.class) plotClusters(ybinpaths, p.class)
3M Markov mixture model for clustering pathways
pathCluster(ybinpaths, M, iter = 1000)
pathCluster(ybinpaths, M, iter = 1000)
ybinpaths |
The training paths computed by |
M |
The number of clusters. |
iter |
The maximum number of EM iterations. |
A list with the following items:
h |
The posterior probabilities that each path belongs to each cluster. |
labels |
The cluster membership labels. |
theta |
The probabilities of each gene for each cluster. |
proportions |
The mixing proportions of each path. |
likelihood |
The likelihood convergence history. |
params |
The specific parameters used. |
Ichigaku Takigawa
Timothy Hancock
Mamitsuka, H., Okuno, Y., and Yamaguchi, A. 2003. Mining biologically active patterns in metabolic pathways using microarray expression profiles. SIGKDD Explor. News l. 5, 2 (Dec. 2003), 113-121.
Other Path clustering & classification methods:
pathClassifier()
,
pathsToBinary()
,
plotClassifierROC()
,
plotClusterMatrix()
,
plotPathClassifier()
,
plotPathCluster()
,
predictPathClassifier()
,
predictPathCluster()
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=8) ## Convert paths to binary matrix. ybinpaths <- pathsToBinary(ranked.p) p.cluster <- pathCluster(ybinpaths, M=2) plotClusters(ybinpaths, p.cluster)
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=8) ## Convert paths to binary matrix. ybinpaths <- pathsToBinary(ranked.p) p.cluster <- pathCluster(ybinpaths, M=2) plotClusters(ybinpaths, p.cluster)
Given a weighted igraph object, path ranking finds a set of node/edge sequences (paths) to
maximize the sum of edge weights.
pathRanker(method="prob.shortest.path")
extracts the K most probable paths within
a weighted network.
pathRanker(method="pvalue")
extracts a list of paths whose sum of edge weights are
significantly higher than random paths of the same length.
pathRanker( graph, method = "prob.shortest.path", start, end, verbose = TRUE, ... )
pathRanker( graph, method = "prob.shortest.path", start, end, verbose = TRUE, ... )
graph |
A weighted igraph object. Weights must be in |
method |
Which path ranking method to use. |
start |
A list of start vertices, given by their vertex id. |
end |
A list of terminal vertices, given by their vertex id. |
verbose |
Whether to display the progress of the function. |
... |
Method-specific parameters. See Details section. |
The input here is graph
. A weight must be assigned to each edge. Bootstrapped Pearson correlation edge weights
can be assigned to each edge by assignEdgeWeights
. However the specification of the edge weight is flexible
with the condition that increasing values indicate stronger relationships between vertices.
pathRanker(method="prob.shortest.path")
finds the K most probable loopless paths given a weighted network.
Before the paths are ranked the edge weights are converted into probabilistic edge weights using the Empirical
Cumulative Distribution (ECDF) over all edge weights. This is called ECDF edge weight. The ECDF edge weight
serves as a probabilistic rank of the most important gene-gene interactions. The probabilistic nature of the ECDF
edge weights allow for a significance test to determine if a path contains any functional structure or is
simply a random walk. The probability of a path is simily the product of all ECDF weights along the path.
This is computed as a sum of the logs of the ECDF edge weights.
The follwing arguments can be passed to pathRanker(method="prob.shortest.path")
:
K
Maximum number of paths to extract. Defaults to 10.
minPathSize
The minimum number of edges for each extracted path. Defualts to 1.
normalize
Specify if you want to normalize the probabilistic edge weights (across different labels) before extracting the paths. Defaults to TRUE.
pathRanker(method="pvalue")
is deprecated. Please use prob.shortest.path
instead.
A list of paths where each path has the following items:
gene |
The ordered sequence of genes visited along the path. |
compounds |
The ordered sequence of compounds visited along the path. |
weights |
The ordered sequence of the log(ECDF edge weights) along the path. |
distance |
The sum of the log(ECDF edge weights) along each path. (a sum of logs is a product) |
Timothy Hancock, Ichigaku Takigawa, Nicolas Wicker and Ahmed Mohamed
getPathsAsEIDs, extractPathNetwork
Other Path ranking methods:
extractPathNetwork()
,
getPathsAsEIDs()
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=6)
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=6)
Converts the result from pathRanker into something suitable for pathClassifier or pathCluster.
pathsToBinary(ypaths)
pathsToBinary(ypaths)
ypaths |
The result of |
Converts a set of pathways from pathRanker
into a list of binary pathway matrices. If the pathways are grouped by a response label then the
pathsToBinary returns a list labeled by response class where each element is the binary
pathway matrix for each class. If the pathways are from pathRanker
then a list wiht
a single element containing the binary pathway matrix is returned. To look up the structure of a
specific binary path in the corresponding ypaths
object simply use matrix index by calling
ypaths[[ybinpaths\$pidx[i,]]]
, where i
is the row in the binary paths object you
wish to reference.
A list with the following elements.
paths |
All paths within ypaths converted to a binary string and concatenated into the one matrix. |
y |
The response variable. |
pidx |
An matrix where each row specifies the location of that path within the |
Timothy Hancock and Ichigaku Takigawa
Other Path clustering & classification methods:
pathClassifier()
,
pathCluster()
,
plotClassifierROC()
,
plotClusterMatrix()
,
plotPathClassifier()
,
plotPathCluster()
,
predictPathClassifier()
,
predictPathCluster()
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=6) ## Convert paths to binary matrix. ybinpaths <- pathsToBinary(ranked.p) p.cluster <- pathCluster(ybinpaths, M=3) plotClusters(ybinpaths, p.cluster, col=c("red", "green", "blue") )
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=6) ## Convert paths to binary matrix. ybinpaths <- pathsToBinary(ranked.p) p.cluster <- pathCluster(ybinpaths, M=3) plotClusters(ybinpaths, p.cluster, col=c("red", "green", "blue") )
This function highlighting ranked paths over different network representations, metabolic, reaction and gene networks. The functions finds equivalent paths across different networks and marks them.
plotAllNetworks( paths, metabolic.net = NULL, reaction.net = NULL, gene.net = NULL, path.clusters = NULL, plot.clusters = TRUE, col.palette = palette(), layout = layout.auto, ... )
plotAllNetworks( paths, metabolic.net = NULL, reaction.net = NULL, gene.net = NULL, path.clusters = NULL, plot.clusters = TRUE, col.palette = palette(), layout = layout.auto, ... )
paths |
The result of |
metabolic.net |
A bipartite metabolic network. |
reaction.net |
A reaction network, resulting from |
gene.net |
A gene network, resulting from |
path.clusters |
The result from |
plot.clusters |
Whether to plot clustering information, as generated by |
col.palette |
A color palette, or a palette generating function (ex: col.palette=rainbow ). |
layout |
Either a graph layout function, or a two-column matrix specifiying vertex coordinates. |
... |
Additional arguments passed to |
Highlights the path list over all provided networks.
Ahmed Mohamed
Other Plotting methods:
colorVertexByAttr()
,
layoutVertexByAttr()
,
plotClassifierROC()
,
plotClusterMatrix()
,
plotCytoscapeGML()
,
plotNetwork()
,
plotPathClassifier()
,
plotPaths()
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=6) plotAllNetworks(ranked.p, metabolic.net = ex_sbml, reaction.net = rgraph, vertex.label = "", vertex.size = 4)
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=6) plotAllNetworks(ranked.p, metabolic.net = ex_sbml, reaction.net = rgraph, vertex.label = "", vertex.size = 4)
Diagnostic plots for pathClassifier
.
plotClassifierROC(mix)
plotClassifierROC(mix)
mix |
The result from |
Diagnostic plots of the result from pathClassifier.
itemTopROC curves for the posterior probabilities (mix\$posterior.probs
)
and for each HME3M component (mix\$h
). This gives information about what response
label each relates to. A ROC curve with an AUC < 0.5
relates to y = 0
.
Conversely ROC curves with AUC > 0.5
relate to y = 1
.
itemBottomThe likelihood convergence history for the HME3M model. If the parameters
alpha
or lambda
are set too large then the likelihood may decrease.
Timothy Hancock and Ichigaku Takigawa
Other Path clustering & classification methods:
pathClassifier()
,
pathCluster()
,
pathsToBinary()
,
plotClusterMatrix()
,
plotPathClassifier()
,
plotPathCluster()
,
predictPathClassifier()
,
predictPathCluster()
Other Plotting methods:
colorVertexByAttr()
,
layoutVertexByAttr()
,
plotAllNetworks()
,
plotClusterMatrix()
,
plotCytoscapeGML()
,
plotNetwork()
,
plotPathClassifier()
,
plotPaths()
Plots the structure of all path clusters
plotClusterMatrix( ybinpaths, clusters, col = rainbow(clusters$params$M), grid = TRUE ) plotClusterProbs(clusters, col = rainbow(clusters$params$M)) plotClusters(ybinpaths, clusters, col, ...)
plotClusterMatrix( ybinpaths, clusters, col = rainbow(clusters$params$M), grid = TRUE ) plotClusterProbs(clusters, col = rainbow(clusters$params$M)) plotClusters(ybinpaths, clusters, col, ...)
ybinpaths |
The training paths computed by |
clusters |
The pathway cluster model trained by |
col |
Colors for each path cluster. |
grid |
A logical, whether to add a |
... |
Extra paramaters passed to |
plotClusterMatrix
plots an image of all paths the training dataset. Rows are the paths and columns
are the genes (features) included within each path. Paths are colored according to cluster membership.
plotClusterProbs
The training set posterior probabilities for each path belonging to a 3M component.
plotClusters
: combines the two plots produced by plotClusterProbs
and plotClusterMatrix
.
Ahmed Mohamed
Other Path clustering & classification methods:
pathClassifier()
,
pathCluster()
,
pathsToBinary()
,
plotClassifierROC()
,
plotPathClassifier()
,
plotPathCluster()
,
predictPathClassifier()
,
predictPathCluster()
Other Plotting methods:
colorVertexByAttr()
,
layoutVertexByAttr()
,
plotAllNetworks()
,
plotClassifierROC()
,
plotCytoscapeGML()
,
plotNetwork()
,
plotPathClassifier()
,
plotPaths()
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=8) ## Convert paths to binary matrix. ybinpaths <- pathsToBinary(ranked.p) p.cluster <- pathCluster(ybinpaths, M=2) plotClusters(ybinpaths, p.cluster, col=c("red", "blue") )
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=8) ## Convert paths to binary matrix. ybinpaths <- pathsToBinary(ranked.p) p.cluster <- pathCluster(ybinpaths, M=2) plotClusters(ybinpaths, p.cluster, col=c("red", "blue") )
plotCytoscape
function has been removed because RCytoscape is no longer prensent in Bioconductor.
Future plans will use RCy3 for Cytoscape plotting, once RCy3 is supported on MacOS and Windows.
plotCytoscapeGML exports the network plot in GML format, that can be later imported into Cytoscape
(using "import network from file" option). This fuction is compatible with all Cytoscape versions.
plotCytoscapeGML( graph, file, layout = layout.auto, vertex.size, vertex.label, vertex.shape, vertex.color, edge.color )
plotCytoscapeGML( graph, file, layout = layout.auto, vertex.size, vertex.label, vertex.shape, vertex.color, edge.color )
graph |
An annotated igraph object. |
file |
Output GML file name to which the network plot is exported. |
layout |
Either a graph layout function, or a two-column matrix specifiying vertex coordinates. |
vertex.size |
Vertex size. If missing, the vertex attribute "size" ( V(g)$size) ) will be used. |
vertex.label |
Vertex labels. If missing, the vertex attribute "label" ( V(g)$label) ) will be used. If missing, vertices are labeled by their name. |
vertex.shape |
Vertex shape in one of igraph shapes. If missing, the vertex attribute "shape" ( V(g)$shape) ) will be used. Shapes are converted from igraph convention to Cytoscape convention. "square","rectangle" and "vrectangle" are converted to "RECT", "csquare" and "crectangle" are converted to "ROUND_RECT", all other shapes are considered "ELLIPSE" |
vertex.color |
A color or a list of colors for vertices. Vetices with multiple colors are not supported. If missing, the vertex attribute "color" ( V(g)$color) ) will be used. |
edge.color |
A color or a list of colors for edges. If missing, the edge attribute "color" ( E(g)$color) ) will be used. |
For plotCytoscapeGML
, results are written to file.
Ahmed Mohamed
Other Plotting methods:
colorVertexByAttr()
,
layoutVertexByAttr()
,
plotAllNetworks()
,
plotClassifierROC()
,
plotClusterMatrix()
,
plotNetwork()
,
plotPathClassifier()
,
plotPaths()
data("ex_sbml") rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) v.layout <- layoutVertexByAttr(rgraph, "compartment") v.color <- colorVertexByAttr(rgraph, "compartment") # Export network plot to GML file plotCytoscapeGML(rgraph, file="example.gml", layout=v.layout, vertex.color=v.color, vertex.size=10)
data("ex_sbml") rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) v.layout <- layoutVertexByAttr(rgraph, "compartment") v.color <- colorVertexByAttr(rgraph, "compartment") # Export network plot to GML file plotCytoscapeGML(rgraph, file="example.gml", layout=v.layout, vertex.color=v.color, vertex.size=10)
This function is a wrapper function for plot.igraph
, with 2 main additions.
1. Add the ability to color vertices by their attributes (see examples), accompanied by an inofrmative
legend. 2. Resize vertex.size, edge.arrow.size, label.cex according to the plot size and the size of the
network.
plotNetwork( graph, vertex.color, col.palette = palette(), layout = layout.auto, legend = TRUE, ... )
plotNetwork( graph, vertex.color, col.palette = palette(), layout = layout.auto, legend = TRUE, ... )
graph |
An annotated igraph object. |
vertex.color |
A list of colors for vertices, or an attribute names (ex: "pathway") by which vertices
will be colored. Complex attributes, where a vertex belongs to more than one group, are supported. This can
also be the output of |
col.palette |
A color palette, or a palette generating function (ex: col.palette=rainbow ). |
layout |
Either a graph layout function, or a two-column matrix specifiying vertex coordinates. |
legend |
Wheter to plot a legend. The legend is only plotted if vertices are colored by attribute values. |
... |
Additional arguments passed to |
Produces a plot of the network.
Ahmed Mohamed
Other Plotting methods:
colorVertexByAttr()
,
layoutVertexByAttr()
,
plotAllNetworks()
,
plotClassifierROC()
,
plotClusterMatrix()
,
plotCytoscapeGML()
,
plotPathClassifier()
,
plotPaths()
data("ex_kgml_sig") plotNetwork(ex_kgml_sig, vertex.color="pathway") plotNetwork(ex_kgml_sig, vertex.color="pathway", col.palette=heat.colors) plotNetwork(ex_kgml_sig, vertex.color="pathway", col.palette=c("red", "green","blue","grey"))
data("ex_kgml_sig") plotNetwork(ex_kgml_sig, vertex.color="pathway") plotNetwork(ex_kgml_sig, vertex.color="pathway", col.palette=heat.colors) plotNetwork(ex_kgml_sig, vertex.color="pathway", col.palette=c("red", "green","blue","grey"))
Plots the structure of specified path found by pathClassifier.
plotPathClassifier(ybinpaths, obj, m, tol = NULL)
plotPathClassifier(ybinpaths, obj, m, tol = NULL)
ybinpaths |
The training paths computed by |
obj |
The pathClassifier |
m |
The path component to view. |
tol |
A tolerance for 3M parameter |
Produces a plot of the paths with the path probabilities and prediction probabilities and ROC curve overlayed.
Center Plot |
An image of all paths the training dataset. Rows are the paths and columns are the genes (vertices) included within each pathway. A colour within image indicates if a particular gene (vertex) is included within a specific path. Colours flag whether a path belongs to the current HME3M component (P > 0.5). |
Center Right |
The training set posterior probabilities for each path belonging to the current 3M component. |
Center Top |
The ROC curve for this HME3M component. |
Top Bar Plots |
|
Timothy Hancock and Ichigaku Takigawa
Other Path clustering & classification methods:
pathClassifier()
,
pathCluster()
,
pathsToBinary()
,
plotClassifierROC()
,
plotClusterMatrix()
,
plotPathCluster()
,
predictPathClassifier()
,
predictPathCluster()
Other Plotting methods:
colorVertexByAttr()
,
layoutVertexByAttr()
,
plotAllNetworks()
,
plotClassifierROC()
,
plotClusterMatrix()
,
plotCytoscapeGML()
,
plotNetwork()
,
plotPaths()
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=6) ## Convert paths to binary matrix. ybinpaths <- pathsToBinary(ranked.p) p.class <- pathClassifier(ybinpaths, target.class = "BCR/ABL", M = 3) ## Plotting the classifier results. plotClassifierROC(p.class) plotClusters(ybinpaths, p.class)
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=6) ## Convert paths to binary matrix. ybinpaths <- pathsToBinary(ranked.p) p.class <- pathClassifier(ybinpaths, target.class = "BCR/ABL", M = 3) ## Plotting the classifier results. plotClassifierROC(p.class) plotClusters(ybinpaths, p.class)
Plots the structure of specified path found by pathCluster.
plotPathCluster(ybinpaths, clusters, m, tol = NULL)
plotPathCluster(ybinpaths, clusters, m, tol = NULL)
ybinpaths |
The training paths computed by |
clusters |
The pathway cluster model trained by |
m |
The path cluster to view. |
tol |
A tolerance for 3M parameter |
Produces a plot of the paths with the path probabilities and cluster membership probabilities.
Center Plot |
An image of all paths the training dataset. Rows are the paths and columns are the genes (features) included within each path. |
Right |
The training set posterior probabilities for each path belonging to the current 3M component. |
Top Bar Plots |
|
Timothy Hancock and Ichigaku Takigawa
Other Path clustering & classification methods:
pathClassifier()
,
pathCluster()
,
pathsToBinary()
,
plotClassifierROC()
,
plotClusterMatrix()
,
plotPathClassifier()
,
predictPathClassifier()
,
predictPathCluster()
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=8) ## Convert paths to binary matrix. ybinpaths <- pathsToBinary(ranked.p) p.cluster <- pathCluster(ybinpaths, M=2) plotPathCluster(ybinpaths, p.cluster, m=2, tol=0.05)
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=8) ## Convert paths to binary matrix. ybinpaths <- pathsToBinary(ranked.p) p.cluster <- pathCluster(ybinpaths, M=2) plotPathCluster(ybinpaths, p.cluster, m=2, tol=0.05)
This function plots a network highlighting ranked paths. If path.clusters
are provided,
paths in the same cluster are assigned similar colors.
plotPaths( paths, graph, path.clusters = NULL, col.palette = palette(), layout = layout.auto, ... )
plotPaths( paths, graph, path.clusters = NULL, col.palette = palette(), layout = layout.auto, ... )
paths |
The result of |
graph |
An annotated igraph object. |
path.clusters |
The result from |
col.palette |
A color palette, or a palette generating function (ex: col.palette=rainbow ). |
layout |
Either a graph layout function, or a two-column matrix specifiying vertex coordinates. |
... |
Additional arguments passed to |
Produces a plot of the network with paths highlighted. If paths are computed for several labels (sample categories), a plot is created for each label.
Ahmed Mohamed
Other Plotting methods:
colorVertexByAttr()
,
layoutVertexByAttr()
,
plotAllNetworks()
,
plotClassifierROC()
,
plotClusterMatrix()
,
plotCytoscapeGML()
,
plotNetwork()
,
plotPathClassifier()
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=6) ## Plot paths. plotPaths(ranked.p, rgraph) ## Convert paths to binary matrix, build a classifier. ybinpaths <- pathsToBinary(ranked.p) p.class <- pathClassifier(ybinpaths, target.class = "BCR/ABL", M = 3) ## Plotting with clusters, on a metabolic graph. plotPaths(ranked.p, ex_sbml, path.clusters=p.class)
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=6) ## Plot paths. plotPaths(ranked.p, rgraph) ## Convert paths to binary matrix, build a classifier. ybinpaths <- pathsToBinary(ranked.p) p.class <- pathClassifier(ybinpaths, target.class = "BCR/ABL", M = 3) ## Plotting with clusters, on a metabolic graph. plotPaths(ranked.p, ex_sbml, path.clusters=p.class)
Predicts new paths given a pathClassifier model.
predictPathClassifier(mix, newdata)
predictPathClassifier(mix, newdata)
mix |
The result from |
newdata |
A data.frame containing the new paths to be classified. |
A list with the following elements.
h |
The posterior probabilities for each HME3M component. |
posterior.probs |
The posterior probabilities for HME3M model to classify the response. |
label |
A vector indicating the HME3M cluster membership. |
component |
The HME3M component membership for each pathway. |
path.probabilities |
The 3M path probabilities. |
plr.probabilities |
The PLR predictions for each component. |
Timothy Hancock and Ichigaku Takigawa
Other Path clustering & classification methods:
pathClassifier()
,
pathCluster()
,
pathsToBinary()
,
plotClassifierROC()
,
plotClusterMatrix()
,
plotPathClassifier()
,
plotPathCluster()
,
predictPathCluster()
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=6) ## Convert paths to binary matrix. ybinpaths <- pathsToBinary(ranked.p) p.class <- pathClassifier(ybinpaths, target.class = "BCR/ABL", M = 3) ## Just an example of how to predict cluster membership pclass.pred <- predictPathCluster(p.class, ybinpaths$paths)
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", y=factor(colnames(ex_microarray)), bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=6) ## Convert paths to binary matrix. ybinpaths <- pathsToBinary(ranked.p) p.class <- pathClassifier(ybinpaths, target.class = "BCR/ABL", M = 3) ## Just an example of how to predict cluster membership pclass.pred <- predictPathCluster(p.class, ybinpaths$paths)
Predicts new paths given a pathCluster model.
predictPathCluster(pfit, newdata)
predictPathCluster(pfit, newdata)
pfit |
The pathway cluster model trained by |
newdata |
The binary pathway dataset to be assigned a cluster label. |
A list with the following elements:
labels |
a vector indicating the 3M cluster membership. |
posterior.probs |
a matrix of posterior probabilities for each path belonging to each cluster. |
Ichigaku Takigawa
Timothy Hancock
Other Path clustering & classification methods:
pathClassifier()
,
pathCluster()
,
pathsToBinary()
,
plotClassifierROC()
,
plotClusterMatrix()
,
plotPathClassifier()
,
plotPathCluster()
,
predictPathClassifier()
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=8) ## Convert paths to binary matrix. ybinpaths <- pathsToBinary(ranked.p) p.cluster <- pathCluster(ybinpaths, M=2) ## just an example of how to predict cluster membership. pclust.pred <- predictPathCluster(p.cluster,ybinpaths$paths)
## Prepare a weighted reaction network. ## Conver a metabolic network to a reaction network. data(ex_sbml) # bipartite metabolic network of Carbohydrate metabolism. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ## Assign edge weights based on Affymetrix attributes and microarray dataset. # Calculate Pearson's correlation. data(ex_microarray) # Part of ALL dataset. rgraph <- assignEdgeWeights(microarray = ex_microarray, graph = rgraph, weight.method = "cor", use.attr="miriam.uniprot", bootstrap = FALSE) ## Get ranked paths using probabilistic shortest paths. ranked.p <- pathRanker(rgraph, method="prob.shortest.path", K=20, minPathSize=8) ## Convert paths to binary matrix. ybinpaths <- pathsToBinary(ranked.p) p.cluster <- pathCluster(ybinpaths, M=2) ## just an example of how to predict cluster membership. pclust.pred <- predictPathCluster(p.cluster,ybinpaths$paths)
Internal method to register memery errors, caused by compiled code. This method is used only by the package, and should not be invoked by users.
registerMemoryErr(method)
registerMemoryErr(method)
method |
The mathod which generated the error. |
Ahmed Mohamed
This function allows users to replace vertex ids with another attribute, calculating connectivities based on the new attribute.
reindexNetwork(graph, v.attr)
reindexNetwork(graph, v.attr)
graph |
An annotated igraph object. |
v.attr |
Name of the attribute to use as vertex ids. |
This functions can be very useful when merging networks constructed from different databases. For example, to match a network created from Reactome to a KEGG network, you can reindex the vertices by "miriam.kegg.compound" attribute. Another usage is to remove duplicated vertices (in case of different subcellular compartments, for example). if a network has ATP_membrane & ATP_cytoplasm vertices, reindexing by chemical name will collapse them into one 'ATP' vertex.
A new graph with vertices expanded.
Ahmed Mohamed
Other Network processing methods:
expandComplexes()
,
makeMetaboliteNetwork()
,
makeReactionNetwork()
,
rmSmallCompounds()
,
simplifyReactionNetwork()
,
vertexDeleteReconnect()
## Make a gene network from a reaction network. data(ex_sbml) # A bipartite metbaolic network. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ggraph <- makeGeneNetwork(rgraph) ## Expand vertices into their contituent genes. data(ex_kgml_sig) # Ras and chemokine signaling pathways in human ggraph <- expandComplexes(ex_kgml_sig, v.attr = "miriam.ncbigene", keep.parent.attr= c("^pathway", "^compartment")) ## Create a separate vertex for each compartment. This is useful in duplicating ## metabolite vertices in a network. ## Not run: graph <- expandComplexes(graph, v.attr = "compartment", keep.parent.attr = "all", expansion.method = "duplicate", missing.method = "keep") ## End(Not run)
## Make a gene network from a reaction network. data(ex_sbml) # A bipartite metbaolic network. rgraph <- makeReactionNetwork(ex_sbml, simplify=TRUE) ggraph <- makeGeneNetwork(rgraph) ## Expand vertices into their contituent genes. data(ex_kgml_sig) # Ras and chemokine signaling pathways in human ggraph <- expandComplexes(ex_kgml_sig, v.attr = "miriam.ncbigene", keep.parent.attr= c("^pathway", "^compartment")) ## Create a separate vertex for each compartment. This is useful in duplicating ## metabolite vertices in a network. ## Not run: graph <- expandComplexes(graph, v.attr = "compartment", keep.parent.attr = "all", expansion.method = "duplicate", missing.method = "keep") ## End(Not run)
This function removes uniquitous compounds (metabolites connected to numerous reactions) from a metabolic network.These compounds are reaction cofactors and currency compounds, such as ATP, CO2, etc. A path through these metabolites may not be bioloigcally meaningful. The defualt small compound list is derived from Reactome, containing keeg.compound, pubchem.compound, ChEBI and CAS identifiers.
rmSmallCompounds( graph, method = c("remove", "duplicate"), small.comp.ls = NPMdefaults("small.comp.ls") )
rmSmallCompounds( graph, method = c("remove", "duplicate"), small.comp.ls = NPMdefaults("small.comp.ls") )
graph |
A metabolic network. |
method |
How to handle small compounds. Either simply delete these vertices |
small.comp.ls |
A list of small compounds to be used. |
A modified graph, with the small compounds removed or duplicated.
Ahmed Mohamed
Other Network processing methods:
expandComplexes()
,
makeMetaboliteNetwork()
,
makeReactionNetwork()
,
reindexNetwork()
,
simplifyReactionNetwork()
,
vertexDeleteReconnect()
data(ex_sbml) sbml.removed <- rmSmallCompounds(ex_sbml, method="remove")
data(ex_sbml) sbml.removed <- rmSmallCompounds(ex_sbml, method="remove")
This function takes SBML files as input, and returns either a metabolic or a signaling network as output.
SBML2igraph( filename, parse.as = c("metabolic", "signaling"), miriam.attr = "all", gene.attr, expand.complexes, verbose = TRUE )
SBML2igraph( filename, parse.as = c("metabolic", "signaling"), miriam.attr = "all", gene.attr, expand.complexes, verbose = TRUE )
filename |
A character vector containing the SBML files to be processed. If a directory path is provided, all *.xml and *.sbml files in it and its subdirectories are included. |
parse.as |
Whether to process file into a metabolic or a signaling network. |
miriam.attr |
A list of annotation attributes to be extracted. If |
gene.attr |
An attribute to distinguish |
expand.complexes |
Split protein complexes into individual gene nodes. Ignored if
|
verbose |
Whether to display the progress of the function. |
Users can specify whether files are processes as metabolic or signaling networks.
Metabolic networks are given as bipartite graphs, where metabolites and reactions represent
vertex types. This is constructed from ListOfReactions
in SBML file, connecting them
to their corresponding substrates and products (ListOfSpecies
). Each reaction vertex has genes
attribute,
listing all modifiers
of this reaction. As a general rule, reactions inherit all annotation
attributes of its catalyzig genes.
Signaling network have genes as vertices and edges represent interactions. Since SBML format may
represent singling events as reaction
, all species are assumed to be genes (rather than small
molecules). For a simple path S0 -> R1 -> S1
, in signaling network, the path will be
S0 -> M(R1) -> S1
where M(R1)
is R1 modifier(s). To ditiguish gene species from small
molecules, user can provide gene.attr
(for example: miriam.uniprot
or miriam.ncbigene
)
where only annotated species are considered genes.
All annotation attributes written according to MIRIAM guidlines (either urn:miriam:xxx:xxx
or
http://identifiers.org/xxx/xxx
) are etxracted by default. Non-conforming attributes can be extracted
by specifying miriam.attr
.
To generate a genome scale network, simply provide a list of files to be parsed, or put all
file in a directory, as pass the directory path as filename
Note: This function requires libSBML installed (Please see the installation instructions in the Vignette). Some SBML level-3 files may requires additional libraries also (An infomative error will be displayed when parsing such files). Please visit http://sbml.org/Documents/Specifications/SBML_Level_3/Packages for more information.
An igraph object, representing a metbolic or a signaling network.
Ahmed Mohamed
Other Database extraction methods:
KGML2igraph()
,
biopax2igraph()
if(is.loaded("readsbmlfile")){ # This is false if libSBML wasn't available at installation. filename <- system.file("extdata", "porphyrin.sbml", package="NetPathMiner") # Process SBML file as a metabolic network g <- SBML2igraph(filename) plotNetwork(g) # Process SBML file as a signaling network g <- SBML2igraph(filename, parse.as="signaling", gene.attr="miriam.uniprot",expand.complexes=TRUE) dev.new() plotNetwork(g) }
if(is.loaded("readsbmlfile")){ # This is false if libSBML wasn't available at installation. filename <- system.file("extdata", "porphyrin.sbml", package="NetPathMiner") # Process SBML file as a metabolic network g <- SBML2igraph(filename) plotNetwork(g) # Process SBML file as a signaling network g <- SBML2igraph(filename, parse.as="signaling", gene.attr="miriam.uniprot",expand.complexes=TRUE) dev.new() plotNetwork(g) }
This function removes reaction vertices with no gene annotations as indicated by the parameter
gene.attr
, and connect their neighbour vertices to preserve graph connectivity. This is
particularly meaningful when reactions are translocation or spontaneous reactions,
which are not catalysed by genes.
simplifyReactionNetwork( reaction.graph, gene.attr = "genes", remove.missing.genes = TRUE, reconnect.threshold = vcount(reaction.graph) )
simplifyReactionNetwork( reaction.graph, gene.attr = "genes", remove.missing.genes = TRUE, reconnect.threshold = vcount(reaction.graph) )
reaction.graph |
A reaction network. |
gene.attr |
The attribute to be considered as "genes". Reactions missing this annotation, will be removed. |
remove.missing.genes |
If |
reconnect.threshold |
An argument passed to |
A simplified reaction network.
Ahmed Mohamed
Other Network processing methods:
expandComplexes()
,
makeMetaboliteNetwork()
,
makeReactionNetwork()
,
reindexNetwork()
,
rmSmallCompounds()
,
vertexDeleteReconnect()
data(ex_sbml) rgraph <- makeReactionNetwork(ex_sbml, simplify=FALSE) ## Removes all reaction nodes with no annotated genes. rgraph <- simplifyReactionNetwork(rgraph, remove.missing.genes=TRUE)
data(ex_sbml) rgraph <- makeReactionNetwork(ex_sbml, simplify=FALSE) ## Removes all reaction nodes with no annotated genes. rgraph <- simplifyReactionNetwork(rgraph, remove.missing.genes=TRUE)
These functions deals with conforming with MIRIAM annotation guidelines, conversion and mapping between MIRIAM identifiers.
stdAttrNames(graph, return.value = c("matches", "graph")) fetchAttribute( graph, organism = "Homo sapiens", target.attr, source.attr, bridge.web = NPMdefaults("bridge.web") )
stdAttrNames(graph, return.value = c("matches", "graph")) fetchAttribute( graph, organism = "Homo sapiens", target.attr, source.attr, bridge.web = NPMdefaults("bridge.web") )
graph |
An annotated igraph object. |
return.value |
Specify whether to return the names of matched standard annotations, or modify the graph attribute names to match the standards. |
organism |
The latin name of the organism (Case-sensitive). |
target.attr |
The target annotation, given as MIRIAM standard in the format |
source.attr |
The source annotation attribute from |
bridge.web |
The base URL for Brigde Database webservices. |
For stdAttrNames
, matches
gives the original attribute names and their MIRIAM version.
Since this is done by simple text matching, mismatches may occur for ambiguous annotations (such as GO, EC number).
graph
returns the input graph with attribute names standardized.
For fetchAttribute
, the input graph
with the fetched attribute mapped to vertices.
Ahmed Mohamed
Other Attribute handling methods:
getAttrStatus()
data(ex_kgml_sig) # Ras and chemokine signaling pathways in human ## Modify attribute names to match MIRIAM standard annotations. graph <- stdAttrNames(ex_kgml_sig, "graph") # Use Attribute fetcher to get affymetrix probeset IDs for network vertices. ## Not run: graph <- fetchAttribute(graph, organism="Homo sapiens", target.attr="miriam.affy.probeset") ## End(Not run)
data(ex_kgml_sig) # Ras and chemokine signaling pathways in human ## Modify attribute names to match MIRIAM standard annotations. graph <- stdAttrNames(ex_kgml_sig, "graph") # Use Attribute fetcher to get affymetrix probeset IDs for network vertices. ## Not run: graph <- fetchAttribute(graph, organism="Homo sapiens", target.attr="miriam.affy.probeset") ## End(Not run)
Converts an annotated igraph object to graphNEL
toGraphNEL(graph, export.attr = "")
toGraphNEL(graph, export.attr = "")
graph |
An annotated igraph object.. |
export.attr |
A |
A graphNEL object.
Ahmed Mohamed
data(ex_kgml_sig) # Ras and chemokine signaling pathways in human graphNEL <- toGraphNEL(ex_kgml_sig, export.attr="^miriam.")
data(ex_kgml_sig) # Ras and chemokine signaling pathways in human graphNEL <- toGraphNEL(ex_kgml_sig, export.attr="^miriam.")
This function removes vertices given as vids
and connects their neighbours as
long as the shortest path beween the neighbours are below the reconnect.threshold
.
vertexDeleteReconnect( graph, vids, reconnect.threshold = vcount(graph), copy.attr = NULL )
vertexDeleteReconnect( graph, vids, reconnect.threshold = vcount(graph), copy.attr = NULL )
graph |
A reaction network. |
vids |
Vertex ids to be removed. |
reconnect.threshold |
If the shortest path between vertices is larger than this threshold, they are not reconnected. |
copy.attr |
A function, or a list of functions, combine edge attributes. Edge attributes of new edges (between reconnected neighbours) are obtained by combining original edges attributes along the shortest path between reconnected neighbors. |
A modified graph.
Ahmed Mohamed
Other Network processing methods:
expandComplexes()
,
makeMetaboliteNetwork()
,
makeReactionNetwork()
,
reindexNetwork()
,
rmSmallCompounds()
,
simplifyReactionNetwork()
## Remove all reaction vertices from a bipartite metabolic network ## keeping only metabolite vertices. data(ex_sbml) graph <- vertexDeleteReconnect(ex_sbml, vids=which(V(ex_sbml)$reactions))
## Remove all reaction vertices from a bipartite metabolic network ## keeping only metabolite vertices. data(ex_sbml) graph <- vertexDeleteReconnect(ex_sbml, vids=which(V(ex_sbml)$reactions))