Title: | GRAPH Interaction from pathway Topological Environment |
---|---|
Description: | Graph objects from pathway topology derived from KEGG, Panther, PathBank, PharmGKB, Reactome SMPDB and WikiPathways databases. |
Authors: | Gabriele Sales [cre], Enrica Calura [aut], Chiara Romualdi [aut] |
Maintainer: | Gabriele Sales <[email protected]> |
License: | AGPL-3 |
Version: | 1.53.0 |
Built: | 2024-12-01 03:44:24 UTC |
Source: | https://github.com/bioc/graphite |
Converts a PathwayList
into a list of Pathway
s.
## S3 method for class 'PathwayList' as.list(x, ...)
## S3 method for class 'PathwayList' as.list(x, ...)
x |
a |
... |
extra arguments to as.list |
A list of pathways.
Gabriele Sales
as.list(pathways("hsapiens", "kegg"))
as.list(pathways("hsapiens", "kegg"))
This function creates a new object of type Pathway
given a data
frame describing its edges.
buildPathway(id, title, species, database, proteinEdges, metaboliteEdges = NULL, mixedEdges = NULL, timestamp = NULL)
buildPathway(id, title, species, database, proteinEdges, metaboliteEdges = NULL, mixedEdges = NULL, timestamp = NULL)
id |
the pathway identifier. |
title |
the title of the pathway. |
species |
the species the pathway belongs to. |
database |
the name of the database the pathway derives from. |
proteinEdges |
a data.frame of edges between proteins (or genes). Must have the following columns: src_type, src, dest_type, dest, direction and type. Direction must be one of the two strings: "directed" or "undirected". |
metaboliteEdges |
interactions between metabolites. Can be |
mixedEdges |
interactions between metabolites and proteins. Can be |
timestamp |
when the pathway was annotated, by default the time
|
A new Pathway
instance.
edges <- data.frame(src_type = "ENTREZID", src="672", dest_type = "ENTREZID", dest="7157", direction="undirected", type="binding") pathway <- buildPathway("#1", "example", "hsapiens", "database", edges) # Example with metabolites: edges <- data.frame(src_type = "ENTREZID", src="672", dest_type = "ENTREZID", dest="7157", direction="undirected", type="binding") mixed <- data.frame(src_type = "CHEBI", src="77750", dest_type = "ENTREZID", dest="7157", direction="undirected", type="binding") pathway <- buildPathway("#1", "example", "hsapiens", "database", edges, mixedEdges = mixed)
edges <- data.frame(src_type = "ENTREZID", src="672", dest_type = "ENTREZID", dest="7157", direction="undirected", type="binding") pathway <- buildPathway("#1", "example", "hsapiens", "database", edges) # Example with metabolites: edges <- data.frame(src_type = "ENTREZID", src="672", dest_type = "ENTREZID", dest="7157", direction="undirected", type="binding") mixed <- data.frame(src_type = "CHEBI", src="77750", dest_type = "ENTREZID", dest="7157", direction="undirected", type="binding") pathway <- buildPathway("#1", "example", "hsapiens", "database", edges, mixedEdges = mixed)
Converts the node identifiers of pathways.
If the option Ncpus
is set to a value larger than 1 and the package
parallel
is installed, the conversion procedure will automatically
use multiple cores.
convertIdentifiers(x, to)
convertIdentifiers(x, to)
x |
can be a list of pathways or a single pathway |
to |
a string describing the type of the identifier.
Can assume the values |
A Pathway
object.
r <- pathways("hsapiens", "reactome") convertIdentifiers(r$`mTORC1-mediated signalling`, "symbol")
r <- pathways("hsapiens", "reactome") convertIdentifiers(r$`mTORC1-mediated signalling`, "symbol")
Renders the topology of a pathway as a Cytoscape graph.
cytoscapePlot(pathway, ..., cy.ver = 3)
cytoscapePlot(pathway, ..., cy.ver = 3)
pathway |
a |
... |
optional arguments forwarded to |
cy.ver |
select a Cytoscape version. Only version 3 is supported in this release. |
Requires the RCy3
package.
An invisible list with two items:
graph |
the |
suid |
the RCy3 network SUID. |
## Not run: r <- pathways() cytoscapePlot(convertIdentifiers(reactome$`Unwinding of DNA`, "symbol")) ## End(Not run)
## Not run: r <- pathways() cytoscapePlot(convertIdentifiers(reactome$`Unwinding of DNA`, "symbol")) ## End(Not run)
"Pathway"
A biological pathway.
A Pathway
instance actually stores multiple variants of the same biological data.
This is the list of included variants:
proteins
: includes only interactions among proteins;
metabolites
: includes only interactions among metabolites;
mixed
: includes all available interactions.
pathwayId(p)
:Returns the native ID of the pathway.
pathwayTitle(p)
:Returns the title of the pathway.
pathwayDatabase(p)
:Returns the name of the database the pathway was derived from.
pathwaySpecies(p)
:Returns the name of the species in which the pathway was annotated.
pathwayTimestamp(p)
:Returns the date of pathway data retrieval.
pathwayURL(p)
:Returns the URL of the pathway in its original database, if available.
convertIdentifiers(p, to)
:Returns a new pathway using a different type of node identifiers.
edges(p, which = c("proteins", "metabolites", "mixed"),
stringsAsFactors = TRUE)
:Returns a data.frame describing the edges of this pathway.
The option which
selects the desired pathway variant (see section "Variants" above).
If stringsAsFactors
is TRUE
, strings are converted to factors.
nodes(p, which = c("proteins", "metabolites", "mixed"))
:Returns the names of the nodes belonging to this pathway.
The option which
selects the desired pathway variant (see section "Variants" above).
plot(p)
:Shows the pathway topology in Cytoscape.
runClipper(p, expr, classes, method, ...)
:Runs a clipper
analysis over the pathway.
runTopologyGSA(p, test, exp1, exp2, alpha, ...)
:Runs a topologyGSA
analysis over the pathway.
Gabriele Sales
reactome <- pathways("hsapiens", "reactome") pathway <- reactome[[1]] pathwayTitle(pathway) pathwaySpecies(pathway) nodes(pathway) edges(pathway)
reactome <- pathways("hsapiens", "reactome") pathway <- reactome[[1]] pathwayTitle(pathway) pathwaySpecies(pathway) nodes(pathway) edges(pathway)
Obtains the list of pathway databases available through graphite
.
pathwayDatabases()
pathwayDatabases()
Returns a data.frame
with two columns: species
and
database
.
Gabriele Sales
pathwayDatabases()
pathwayDatabases()
Builds a graphNEL
object representing the topology of a pathway.
pathwayGraph(pathway, which = "proteins", edge.types = NULL)
pathwayGraph(pathway, which = "proteins", edge.types = NULL)
pathway |
a |
which |
the pathway variant you want. See |
edge.types |
keep only the edges maching the type names in this vector. |
A graphNEL
object.
r <- pathways("hsapiens", "reactome") pathwayGraph(r$`mTORC1-mediated signalling`, edge.types="Binding")
r <- pathways("hsapiens", "reactome") pathwayGraph(r$`mTORC1-mediated signalling`, edge.types="Binding")
"PathwayList"
A collection of pathways from a single database.
Class "Pathways"
, directly.
l[i]
returns a selection of the pathways contained in the pathway list.
l[[i]]
gives access to one of the pathways contained in the pathway list.
l$`title`
loads a pathways by its title.
convertIdentifiers(l, to)
returns a new list of pathways using a different type of node identifiers.
length(l)
returns the number of pathways contained in the list.
names(l)
returns the titles of the pathways contained in the list.
prepareSPIA(l, pathwaySetName, print.names=FALSE)
prepares the pathways for a SPIA analysis.
runClipper(l, expr, classes, method, maxNodes=150, ...)
runs a clipper
analysis over all the pathways in the list.
runTopologyGSA(l, test, exp1, exp2, alpha, maxNodes=150, ...)
runs a topologyGSA
analysis over all the pathways in the list.
Gabriele Sales
Retrieve a list of pathways from a database for a given species.
graphite currently supports the following databases:
Call the pathwayDatabase
function for more details.
pathways(species, database)
pathways(species, database)
species |
one of the supported species |
database |
the name of the pathway database |
A PathwayList
object.
pathways("hsapiens", "reactome")
pathways("hsapiens", "reactome")
"Pathways"
A virtual class acting as a common parent to all other classes representing pathway databases.
A virtual Class: No objects may be created from it.
No methods defined with class "Pathways" in the signature.
Gabriele Sales
Prepare pathway dataset needed by runSPIA. See runSPIA
and
spia
for more details.
prepareSPIA(db, pathwaySetName, print.names = FALSE)
prepareSPIA(db, pathwaySetName, print.names = FALSE)
db |
a |
pathwaySetName |
name of the output pathway set. |
print.names |
print pathway names as the conversion advances. |
This function has no return value.
Tarca AL, Draghici S, Khatri P, Hassan SS, Mittal P, Kim JS, Kim CJ, Kusanovic JP, Romero R. A novel signaling pathway impact analysis. Bioinformatics. 2009 Jan 1;25(1):75-82.
Adi L. Tarca, Sorin Draghici, Purvesh Khatri, et. al, A Signaling Pathway Impact Analysis for Microarray Experiments, 2008, Bioinformatics, 2009, 25(1):75-82.
Draghici, S., Khatri, P., Tarca, A.L., Amin, K., Done, A., Voichita, C., Georgescu, C., Romero, R.: A systems biology approach for pathway level analysis. Genome Research, 17, 2007.
Run a topological analysis on an expression dataset using SPIA.
runSPIA(de, all, pathwaySetName, ...)
runSPIA(de, all, pathwaySetName, ...)
de |
A named vector containing log2 fold-changes of the differentially expressed genes. The names of this numeric vector are Entrez gene IDs. |
all |
A vector with the Entrez IDs in the reference set. If the data was obtained from a microarray experiment, this set will contain all genes present on the specific array used for the experiment. This vector should contain all names of the 'de' argument. |
pathwaySetName |
The name of a pathway set created with |
... |
Additional options to pass to |
The spia option "organism" is internally used. It is an error use it in the additional options.
The same of spia, without KEGG links. A data frame containing the ranked pathways and various statistics: pSize is the number of genes on the pathway; NDE is the number of DE genes per pathway; tA is the observed total preturbation accumulation in the pathway; pNDE is the probability to observe at least NDE genes on the pathway using a hypergeometric model; pPERT is the probability to observe a total accumulation more extreme than tA only by chance; pG is the p-value obtained by combining pNDE and pPERT; pGFdr and pGFWER are the False Discovery Rate and respectively Bonferroni adjusted global p-values; and the Status gives the direction in which the pathway is perturbed (activated or inhibited).
Tarca AL, Draghici S, Khatri P, Hassan SS, Mittal P, Kim JS, Kim CJ, Kusanovic JP, Romero R. A novel signaling pathway impact analysis. Bioinformatics. 2009 Jan 1;25(1):75-82.
Adi L. Tarca, Sorin Draghici, Purvesh Khatri, et. al, A Signaling Pathway Impact Analysis for Microarray Experiments, 2008, Bioinformatics, 2009, 25(1):75-82.
Draghici, S., Khatri, P., Tarca, A.L., Amin, K., Done, A., Voichita, C., Georgescu, C., Romero, R.: A systems biology approach for pathway level analysis. Genome Research, 17, 2007.
if (require(SPIA) && require(hgu133plus2.db)) { data(colorectalcancer) top$ENTREZ <- mapIds(hgu133plus2.db, top$ID, "ENTREZID", "PROBEID", multiVals = "first") top <- top[!is.na(top$ENTREZ) & !duplicated(top$ENTREZ), ] top$ENTREZ <- paste("ENTREZID", top$ENTREZ, sep = ":") tg1 <- top[top$adj.P.Val < 0.05, ] DE_Colorectal = tg1$logFC names(DE_Colorectal) <- tg1$ENTREZ ALL_Colorectal <- top$ENTREZ kegg <- pathways("hsapiens", "kegg")[1:20] kegg <- convertIdentifiers(kegg, "ENTREZID") prepareSPIA(kegg, "keggEx") runSPIA(de = DE_Colorectal, all = ALL_Colorectal, "keggEx") unlink("keggExSPIA.RData") }
if (require(SPIA) && require(hgu133plus2.db)) { data(colorectalcancer) top$ENTREZ <- mapIds(hgu133plus2.db, top$ID, "ENTREZID", "PROBEID", multiVals = "first") top <- top[!is.na(top$ENTREZ) & !duplicated(top$ENTREZ), ] top$ENTREZ <- paste("ENTREZID", top$ENTREZ, sep = ":") tg1 <- top[top$adj.P.Val < 0.05, ] DE_Colorectal = tg1$logFC names(DE_Colorectal) <- tg1$ENTREZ ALL_Colorectal <- top$ENTREZ kegg <- pathways("hsapiens", "kegg")[1:20] kegg <- convertIdentifiers(kegg, "ENTREZID") prepareSPIA(kegg, "keggEx") runSPIA(de = DE_Colorectal, all = ALL_Colorectal, "keggEx") unlink("keggExSPIA.RData") }
Use graphical models to test the pathway components highlighting those involved in its deregulation.
If the option Ncpus
is set to a value larger than 1 and the package
parallel
is installed, the conversion procedure will automatically
use multiple cores.
runTopologyGSA(x, test, exp1, exp2, alpha, ...)
runTopologyGSA(x, test, exp1, exp2, alpha, ...)
x |
a |
test |
Either |
exp1 |
Experiment matrix of the first class, genes in columns. |
exp2 |
Experiment matrix of the second class, genes in columns. |
alpha |
Significance level of the test. |
... |
Additional parameters forwarded to When invoked on a |
This function produces a warning and returns NULL when the number of genes in common between the expression matrices and the pathway is less than 3.
See documentation of
pathway.var.test
and
pathway.mean.test
.
Massa MS, Chiogna M, Romualdi C. Gene set analysis exploiting the topology of a pathway. BMC System Biol. 2010 Sep 1;4:121.
if (require(topologyGSA)) { data(examples) colnames(y1) <- paste("SYMBOL", colnames(y1), sep = ":") colnames(y2) <- paste("SYMBOL", colnames(y2), sep = ":") k <- pathways("hsapiens", "kegg") p <- convertIdentifiers(k[["Fc epsilon RI signaling pathway"]], "SYMBOL") runTopologyGSA(p, "var", y1, y2, 0.05) }
if (require(topologyGSA)) { data(examples) colnames(y1) <- paste("SYMBOL", colnames(y1), sep = ":") colnames(y2) <- paste("SYMBOL", colnames(y2), sep = ":") k <- pathways("hsapiens", "kegg") p <- convertIdentifiers(k[["Fc epsilon RI signaling pathway"]], "SYMBOL") runTopologyGSA(p, "var", y1, y2, 0.05) }