Title: | Prediction of pathway-related protein-protein interaction networks |
---|---|
Description: | Package to predict protein-protein interaction (PPI) networks in target organisms for which only a view information about PPIs is available. Path2PPI predicts PPI networks based on sets of proteins which can belong to a certain pathway from well-established model organisms. It helps to combine and transfer information of a certain pathway or biological process from several reference organisms to one target organism. Path2PPI only depends on the sequence similarity of the involved proteins. |
Authors: | Oliver Philipp [aut, cre], Ina Koch [ctb] |
Maintainer: | Oliver Philipp <[email protected]> |
License: | GPL (>= 2) |
Version: | 1.37.0 |
Built: | 2024-10-30 09:19:27 UTC |
Source: | https://github.com/bioc/Path2PPI |
Adds reference species to an object from the class Path2PPI
.
addReference(path2ppi, taxName, taxId, proteins, irefindex, homologs)
addReference(path2ppi, taxName, taxId, proteins, irefindex, homologs)
path2ppi |
An object of the class |
taxName |
A character string giving the taxonomy name. |
taxId |
A character string giving the taxonomy identifier. |
proteins |
Either a character vector with the identifiers of the proteins which are involved in the corresponding pathway or a character vector with the protein names or aliases, respectively, named by the protein identifiers. |
irefindex |
Either a data frame, representing the iRefIndex table of the current reference
species, e.g. loaded previously via |
homologs |
Either a data frame representing the results of the BLAST search (e.g. loaded
previously via |
This method searches for all relevant interactions in the data frame or file
defined in iRefIndex. There are different and often ambiguous protein
identifiers defined in an iRefIndex file, and the putative "major" identifiers
are not necessarily those defined in the corresponding "major" columns "uidA"
and "uidB". Furthermore, iRefIndex also contains protein complexes. Hence,
Path2PPI
applies an advanced search algorithm to automatically find
relevant interactions associated with the pathway or the proteins of interest,
respectively. The user does not have to predefine the identifiers types
(Uniprot, Swissprot, Ensemble etc.), since these types are often unambiguously
assigned. The algorithm searches for each identifier in 10 columns where any
type of identifier or accession number is defined ("uidA", "altA",
"OriginalReferenceA", "FinalReferenceA", "aliasA", "uidB", "altB",
"OriginalReferenceB", "FinalReferenceB" and "aliasB"). Additionally, it
searches for each complex which contains one or more of the predefined
proteins. Subsequently, each homologous relationship which is not relevant for
the previously found interactions is declined. The results of these searches
are centralized in the Path2PPI
object and can be visualized using the
appropriate methods (e.g. showReferences
)
An object from the class Path2PPI
with attached reference species.
Oliver Philipp [email protected]
showReferences
, removeReference
data(ai) #Load test data set ppi <- Path2PPI("Autophagy induction", "Podospora anserina", "5145") ppi ppi <- addReference(ppi, "Homo sapiens", "9606", human.ai.proteins, human.ai.irefindex, pa2human.ai.homologs) ppi <- addReference(ppi, "Saccharomyces cerevisiae (S288c)", "559292", yeast.ai.proteins, yeast.ai.irefindex, pa2yeast.ai.homologs) ppi
data(ai) #Load test data set ppi <- Path2PPI("Autophagy induction", "Podospora anserina", "5145") ppi ppi <- addReference(ppi, "Homo sapiens", "9606", human.ai.proteins, human.ai.irefindex, pa2human.ai.homologs) ppi <- addReference(ppi, "Saccharomyces cerevisiae (S288c)", "559292", yeast.ai.proteins, yeast.ai.irefindex, pa2yeast.ai.homologs) ppi
This data set consists of all data files necessary to predict the putative interactions of the induction step of autophagy in Podospora anserina by means of the corresponding PPIs in human and yeast.
data("ai")
data("ai")
human.ai.irefindex: Data frames with 1694 observations of 54 variables. yeast.ai.irefindex: Data frames with 3840 observations of 54 variables. pa2human.ai.homologs: Data frames with 261 observations of 12 variables. pa2yeast.ai.homologs: Data frames with 98 observations of 12 variables. human.ai.proteins: Named character vector with 5 elements. yeast.ai.proteins: Named character vector with 7 elements.
Data frames human.ai.irefindex
and yeast.ai.irefindex
consists of all relevant interactions of the corresponding iRefIndex
files. The two data frames pa2human.ai.homologs
and
pa2yeast.ai.homologs
are the necessary parts of the result files
from the BLAST searches of the P. anserina proteom against the
proteoms of human and yeast. The named character vectors
human.ai.proteins
and yeast.ai.proteins
consists of the
proteins involved in the induction process of autophagy in human and
yeast.
Four data frames and two named character vectors (see above).
Camacho, C. et al. (2009). BLAST+: architecture and applications. BMC Bioinformatics, 10(1), 421.
Razick, S. et al. (2008). iRefIndex: a consolidated protein interaction database with provenance. BMC Bioinformatics, 9(1), 405.
data(ai)
data(ai)
Get the hybrid network of the previously predicted PPI. The hybrid network consists of all relevant interactions from the reference species, the predicted interactions in the target species and all relevant homologous relationships.
getHybridNetwork(path2ppi, igraph = FALSE)
getHybridNetwork(path2ppi, igraph = FALSE)
path2ppi |
An object of the class |
igraph |
Logical; if |
See igraph
argument.
Oliver Philipp [email protected]
data(ai) #Load test data set ppi <- Path2PPI("Autophagy induction", "Podospora anserina", "5145") ppi <- addReference(ppi, "Homo sapiens", "9606", human.ai.proteins, human.ai.irefindex, pa2human.ai.homologs) ppi <- addReference(ppi, "Saccharomyces cerevisiae (S288c)", "559292", yeast.ai.proteins, yeast.ai.irefindex, pa2yeast.ai.homologs) ppi <- predictPPI(ppi) #Return the hybrid network as data frame hybrid <- getHybridNetwork(ppi) #Return the hybrid network as igraph object hybrid <- getHybridNetwork(ppi,igraph=TRUE)
data(ai) #Load test data set ppi <- Path2PPI("Autophagy induction", "Podospora anserina", "5145") ppi <- addReference(ppi, "Homo sapiens", "9606", human.ai.proteins, human.ai.irefindex, pa2human.ai.homologs) ppi <- addReference(ppi, "Saccharomyces cerevisiae (S288c)", "559292", yeast.ai.proteins, yeast.ai.irefindex, pa2yeast.ai.homologs) ppi <- predictPPI(ppi) #Return the hybrid network as data frame hybrid <- getHybridNetwork(ppi) #Return the hybrid network as igraph object hybrid <- getHybridNetwork(ppi,igraph=TRUE)
Get the predicted PPI of an Path2PPI
object consisting of each
predicted interaction and protein in the target species.
getPPI(path2ppi, raw=FALSE, igraph=FALSE)
getPPI(path2ppi, raw=FALSE, igraph=FALSE)
path2ppi |
An object of the class |
raw |
Logical; if |
igraph |
Logical; if |
See igraph
argument.
Oliver Philipp [email protected]
data(ai) #Load test data set ppi <- Path2PPI("Autophagy induction", "Podospora anserina", "5145") ppi <- addReference(ppi, "Homo sapiens", "9606", human.ai.proteins, human.ai.irefindex, pa2human.ai.homologs) ppi <- addReference(ppi, "Saccharomyces cerevisiae (S288c)", "559292", yeast.ai.proteins, yeast.ai.irefindex, pa2yeast.ai.homologs) ppi <- predictPPI(ppi) #Get the predicted PPI as data frame. network <- getPPI(ppi) #Get the detailed predicted PPI as data frame. network.raw <- getPPI(ppi,raw=TRUE)
data(ai) #Load test data set ppi <- Path2PPI("Autophagy induction", "Podospora anserina", "5145") ppi <- addReference(ppi, "Homo sapiens", "9606", human.ai.proteins, human.ai.irefindex, pa2human.ai.homologs) ppi <- addReference(ppi, "Saccharomyces cerevisiae (S288c)", "559292", yeast.ai.proteins, yeast.ai.irefindex, pa2yeast.ai.homologs) ppi <- predictPPI(ppi) #Get the predicted PPI as data frame. network <- getPPI(ppi) #Get the detailed predicted PPI as data frame. network.raw <- getPPI(ppi,raw=TRUE)
Computes the homology scores based on the BLAST E-value. This function is used
by the predictPPI
method to compute homology scores to decide whether
an interaction in a reference species is adopted to the target species
(see package vignette for a detailed description). It can be used to test
which E-values lead to which scores given a predefined E-value range.
homologyScore(e.value, h.range)
homologyScore(e.value, h.range)
e.value |
One BLAST E-value or a numeric vector with different BLAST E-values |
h.range |
Numeric vector consisting of two values. The first value indicates the lower bound (smallest E-value). Each E-value which is equal or less than this bound is scored with 1. The second value indicates the upper bound (biggest E-value). Each E-value which is equal or greater than this bound is scored with 0. |
Uses a linear function to map the E-value to the range
where
is the lower and
the upper bound:
Numeric vector containing the scores.
Oliver Philipp [email protected]
l <- 1e-100 #lower bound u <- 1e-20 #upper bound h.range <- c(l,u) #define range e.values <- c(1e-20,1e-40,1e-60,1e-80,1e-100) #some BLAST E-values homologyScore(e.values,h.range)
l <- 1e-100 #lower bound u <- 1e-20 #upper bound h.range <- c(l,u) #define range e.values <- c(1e-20,1e-40,1e-60,1e-80,1e-100) #some BLAST E-values homologyScore(e.values,h.range)
"Path2PPI"
An instance of the class Path2PPI
is the major object in the Path2PPI
package. It manages all reference species and the target species.
The prediction algorithm is implemented in this class as well.
Path2PPI(...)
Path2PPI(...)
... |
Argument list (see Note below). |
An instance of the class Path2PPI
.
pathway
:Object of class "character"
targetSpecies
:Object of class ".TargetSpecies"
referenceContainer
:Object of class
".ReferenceContainer"
h.thresh
:Object of class "numeric"
h.range
:Object of class "numeric"
i.thresh
:Object of class "numeric"
consider.complexes
:Object of class "logical"
max.complex.size
:Object of class "numeric"
raw.ppi
:Object of class "data.frame"
ppi
:Object of class "data.frame"
addReference
signature(path2ppi = "Path2PPI")
getHybridNetwork
signature(path2ppi = "Path2PPI")
getPPI
signature(path2ppi = "Path2PPI")
initialize
signature(.Object = "Path2PPI")
plot.Path2PPI
signature(x = "Path2PPI")
predictPPI
signature(path2ppi = "Path2PPI")
removeReference
signature(path2ppi = "Path2PPI")
show
signature(object = "Path2PPI")
showInteraction
signature(path2ppi = "Path2PPI")
showReferences
signature(path2ppi = "Path2PPI")
Arguments to Path2PPI()
and the new
method are obligatory
and must be named if they differ from this order:
pathway | A character string with the name of the pathway which has to be predicted. |
targetName | A character string giving the taxonomy name of the target species. |
targetId | A character string giving the taxonomy identifier of the target species. |
Oliver Philipp [email protected]
ppi <- Path2PPI("Autophagy induction", "Podospora anserina", "5145") ppi
ppi <- Path2PPI("Autophagy induction", "Podospora anserina", "5145") ppi
Plots the predicted PPI in three different ways. Depending on the type
argument it manages the specific layout settings and finally uses the plot
function of the igraph
package.
## S3 method for class 'Path2PPI' plot(x, type = "ppi", multiple.edges = FALSE, scores = FALSE, species.colors = c(), vertices.opacity=0.8, use.identifiers=FALSE, protein.labels = NA, show.legend = TRUE, vertices.coordinates = NA, return.coordinates = FALSE, tkplot=FALSE,...)
## S3 method for class 'Path2PPI' plot(x, type = "ppi", multiple.edges = FALSE, scores = FALSE, species.colors = c(), vertices.opacity=0.8, use.identifiers=FALSE, protein.labels = NA, show.legend = TRUE, vertices.coordinates = NA, return.coordinates = FALSE, tkplot=FALSE,...)
x |
An object from the class |
type |
Character string. Which graph type to plot. "ppi": plots only the predicted PPI. "hybrid": plots the hybrid network which consists of all relevant interactions from the reference species, the predicted interactions in the target species and all relevant homologous relationships. |
multiple.edges |
Logical. Is only considered if |
scores |
Logical. If |
species.colors |
Named vector, to specify the species colors. If no value is given then default colors are used. |
vertices.opacity |
Numeric value between 0 and 1 defining the opacity of the vertices. |
use.identifiers |
Logical. If |
protein.labels |
Named vector to define the labels of the vertices. If no value is given then the protein identifiers are used. The vector does not have to be complete, i.e. not each protein has to be defined. |
show.legend |
Logical. If |
vertices.coordinates |
Data frame containing the coordinates of the vertices. If no value is given
then coordinates are computed using the |
return.coordinates |
Logical. If |
tkplot |
Logical. If |
... |
Additional plotting parameters. |
The argument return.coordinates
only works correctly if
tkplot=FALSE
. If you want to get the coordinates of the tkplot device
use tkplot.getcoords
.
If return.coordinates=TRUE
the coordinates of the vertices are
returned.
If you want to export the plotted graph to postscript you have to consider
that the default font family is set to sans for vertex and edge labels.
Please change the default font family of postscript to sans before you
call the plot method: ps.options(family="sans")
.
Additionally, you have to consider that the default value for
vertices.opacity
is set to 0.8 in order to enhance the
visibility of the graph, since some edges may be hidden by the vertices.
Postscript does not support semi-transparencies. Hence, please change the
vertices.opacity
argument to 1 if you want to export the graph
using postscript.
Oliver Philipp [email protected]
predictPPI
, igraph
for other plotting parameters
data(ai) #Load test data set ppi <- Path2PPI("Autophagy induction", "Podospora anserina", "5145") ppi <- addReference(ppi, "Homo sapiens", "9606", human.ai.proteins, human.ai.irefindex, pa2human.ai.homologs) ppi <- addReference(ppi, "Saccharomyces cerevisiae (S288c)", "559292", yeast.ai.proteins, yeast.ai.irefindex, pa2yeast.ai.homologs) ppi <- predictPPI(ppi,h.range=c(1e-60,1e-20)) #Plot the predicted PPI with the default settings and return #the coordinates of the vertices set.seed(12) coordinates <- plot(ppi, return.coordinates=TRUE) #Plot the predicted PPI and show each underlying reference interaction. #Use different species specific colors. To compare both graphs, #use the coordinates computed before plot(ppi,multiple.edges=TRUE,vertices.coordinates=coordinates) #Plot the corresponding hybrid network with predefined species colors. #Also define some labels for the proteins of the target species. #Keep in mind: You can not use the data in "coordinates" since #the hybrid network consists of more vertices than the default PPI set.seed(40) target.labels<-c("B2AE79"="PaTOR","B2AXK6"="PaATG1", "B2AUW3"="PaATG17","B2AM44"="PaATG11", "B2AQV0"="PaATG13","B2B5M3"="PaVAC8") species.colors <- c("5145"="red","9606"="blue","559292"="green") plot(ppi,type="hybrid",species.colors=species.colors, protein.labels=target.labels)
data(ai) #Load test data set ppi <- Path2PPI("Autophagy induction", "Podospora anserina", "5145") ppi <- addReference(ppi, "Homo sapiens", "9606", human.ai.proteins, human.ai.irefindex, pa2human.ai.homologs) ppi <- addReference(ppi, "Saccharomyces cerevisiae (S288c)", "559292", yeast.ai.proteins, yeast.ai.irefindex, pa2yeast.ai.homologs) ppi <- predictPPI(ppi,h.range=c(1e-60,1e-20)) #Plot the predicted PPI with the default settings and return #the coordinates of the vertices set.seed(12) coordinates <- plot(ppi, return.coordinates=TRUE) #Plot the predicted PPI and show each underlying reference interaction. #Use different species specific colors. To compare both graphs, #use the coordinates computed before plot(ppi,multiple.edges=TRUE,vertices.coordinates=coordinates) #Plot the corresponding hybrid network with predefined species colors. #Also define some labels for the proteins of the target species. #Keep in mind: You can not use the data in "coordinates" since #the hybrid network consists of more vertices than the default PPI set.seed(40) target.labels<-c("B2AE79"="PaTOR","B2AXK6"="PaATG1", "B2AUW3"="PaATG17","B2AM44"="PaATG11", "B2AQV0"="PaATG13","B2B5M3"="PaVAC8") species.colors <- c("5145"="red","9606"="blue","559292"="green") plot(ppi,type="hybrid",species.colors=species.colors, protein.labels=target.labels)
Major method of the Path2PPI
class to predict the final PPI in the
target species using the information available from the stored reference
species. Different values for the arguments of this method can lead to
different PPI networks, differing in the degree of reliability and strictness.
predictPPI(path2ppi, mode="both", h.thresh=1e-05, h.range=c(1e-100, 1e-20), i.thresh=0.7, consider.complexes=FALSE, max.complex.size=5, decline.self.interaction.ref=FALSE, decline.self.interaction.tar=TRUE, verbose=TRUE)
predictPPI(path2ppi, mode="both", h.thresh=1e-05, h.range=c(1e-100, 1e-20), i.thresh=0.7, consider.complexes=FALSE, max.complex.size=5, decline.self.interaction.ref=FALSE, decline.self.interaction.tar=TRUE, verbose=TRUE)
path2ppi |
An object of the class |
mode |
Which interaction from the reference species should be taken into account. "both": both interactors of an interaction has to be in the initial protein list previously inserted by the user (recommended if it is a large network or many proteins were initially defined, respectively). "one": only one of the interactors of each reference interaction has to be in the initial protein list (may lead to very large networks). |
h.thresh |
E-value cutoff at which each homologous relationship definitely will be
declined (see also |
h.range |
Numeric vector consisting of two values. The first value indicates the lower border (smallest E-value). Each E-value which is equal or less than this border is scored with 1 (best). The second value indicates the upper border (biggest E-value). Each E-value which is equal or greater than this border is scored with 0 (worst). |
i.thresh |
Numeric. Threshold for accepted interactions. If the computed prediction score
for an interaction is less than |
consider.complexes |
Logical. If |
max.complex.size |
Numeric. Is only considered if |
decline.self.interaction.ref |
Logical. If |
decline.self.interaction.tar |
Logical. If |
verbose |
Logical. |
Difference of h.thresh
and h.range
: If only one protein in the
target species was found to be homologous to a current reference species
protein and this homology was rated with an E-value which is equal or smaller
than h.thresh
it is scored with 1 (even if the E-value is larger than
the upper border of h.range
). See package vignette for more details.
Use the complex arguments with care, since each complex may lead to a vast amount of interactions, i.e., each protein is considered to interact with each other of this complex; e.g. if there are 10 proteins involved in one complex, this would lead to 10 over 2 = 45 interactions.
An object of the class Path2PPI
with predicted PPI.
Oliver Philipp [email protected]
data(ai) #Load test data set ppi <- Path2PPI("Autophagy induction", "Podospora anserina", "5145") ppi <- addReference(ppi, "Homo sapiens", "9606", human.ai.proteins, human.ai.irefindex, pa2human.ai.homologs) ppi <- addReference(ppi, "Saccharomyces cerevisiae (S288c)", "559292", yeast.ai.proteins, yeast.ai.irefindex, pa2yeast.ai.homologs) #Using the default settings leads to 8 predicted interactions in the #target species ppi <- predictPPI(ppi) #Consider complexes where each complex is allowed to be up to 10 proteins #large. For this smaller pathway only one more interaction was predicted when #considering larger complexes. ppi <- predictPPI(ppi,consider.complexes=TRUE,max.complex.size=10) #We can be less strict and decrease h.range what obviously increases the #number of predicted interactions to 13 ppi <- predictPPI(ppi,h.range=c(1e-60,1e-20))
data(ai) #Load test data set ppi <- Path2PPI("Autophagy induction", "Podospora anserina", "5145") ppi <- addReference(ppi, "Homo sapiens", "9606", human.ai.proteins, human.ai.irefindex, pa2human.ai.homologs) ppi <- addReference(ppi, "Saccharomyces cerevisiae (S288c)", "559292", yeast.ai.proteins, yeast.ai.irefindex, pa2yeast.ai.homologs) #Using the default settings leads to 8 predicted interactions in the #target species ppi <- predictPPI(ppi) #Consider complexes where each complex is allowed to be up to 10 proteins #large. For this smaller pathway only one more interaction was predicted when #considering larger complexes. ppi <- predictPPI(ppi,consider.complexes=TRUE,max.complex.size=10) #We can be less strict and decrease h.range what obviously increases the #number of predicted interactions to 13 ppi <- predictPPI(ppi,h.range=c(1e-60,1e-20))
Remove reference species previously attached to an object from the class
Path2PPI
.
removeReference(path2ppi, species)
removeReference(path2ppi, species)
path2ppi |
An object from the class |
species |
Either a number between 1 and the number of stored reference species or a character string with the taxonomy id of the reference species to remove. |
An object of the class Path2PPI
with removed reference species
species
.
Oliver Philipp [email protected]
data(ai) #Load test data set ppi <- Path2PPI("Autophagy induction", "Podospora anserina", "5145") ppi <- addReference(ppi, "Homo sapiens", "9606", human.ai.proteins, human.ai.irefindex, pa2human.ai.homologs) ppi <- addReference(ppi, "Saccharomyces cerevisiae (S288c)", "559292", yeast.ai.proteins, yeast.ai.irefindex, pa2yeast.ai.homologs) #Remove second reference species ppi <- removeReference(ppi,2) #Remove reference species with taxonomy id "9606" ppi <- removeReference(ppi,"9606")
data(ai) #Load test data set ppi <- Path2PPI("Autophagy induction", "Podospora anserina", "5145") ppi <- addReference(ppi, "Homo sapiens", "9606", human.ai.proteins, human.ai.irefindex, pa2human.ai.homologs) ppi <- addReference(ppi, "Saccharomyces cerevisiae (S288c)", "559292", yeast.ai.proteins, yeast.ai.irefindex, pa2yeast.ai.homologs) #Remove second reference species ppi <- removeReference(ppi,2) #Remove reference species with taxonomy id "9606" ppi <- removeReference(ppi,"9606")
Use showInteraction
to get detailed information about one
interaction of the predicted PPI.
showInteraction(path2ppi, interaction, mode="default", verbose=TRUE)
showInteraction(path2ppi, interaction, mode="default", verbose=TRUE)
path2ppi |
An object from the class |
interaction |
Character vector consisting of the identifiers of the two interactors. |
mode |
Character string. Which information of this interaction is requested. "default": only the predicted interaction and some major information are provided. "detailed": all interactions deduced from each reference species with this interaction is provided. "references": each reference interaction of the current interaction with some major information. "references.detailed": each reference interaction of the current interaction with all available information (extracted from the corresponding iRefIndex data set). |
verbose |
Logical. |
Data frame with the requested information defined in mode
.
Oliver Philipp [email protected]
data(ai) #Load test data set ppi <- Path2PPI("Autophagy induction", "Podospora anserina", "5145") ppi <- addReference(ppi, "Homo sapiens", "9606", human.ai.proteins, human.ai.irefindex, pa2human.ai.homologs) ppi <- addReference(ppi, "Saccharomyces cerevisiae (S288c)", "559292", yeast.ai.proteins, yeast.ai.irefindex, pa2yeast.ai.homologs) ppi <- predictPPI(ppi,h.range=c(1e-60,1e-20)) interaction <- showInteraction(ppi,interaction=c("B2AT71","B2AE79"), mode="detailed") interaction
data(ai) #Load test data set ppi <- Path2PPI("Autophagy induction", "Podospora anserina", "5145") ppi <- addReference(ppi, "Homo sapiens", "9606", human.ai.proteins, human.ai.irefindex, pa2human.ai.homologs) ppi <- addReference(ppi, "Saccharomyces cerevisiae (S288c)", "559292", yeast.ai.proteins, yeast.ai.irefindex, pa2yeast.ai.homologs) ppi <- predictPPI(ppi,h.range=c(1e-60,1e-20)) interaction <- showInteraction(ppi,interaction=c("B2AT71","B2AE79"), mode="detailed") interaction
Get information about the currently stored reference species. If indicated by
returnValue
a data frame - containing information about each protein or
interaction - is provided as well.
showReferences(path2ppi, species = NA, returnValue = NA)
showReferences(path2ppi, species = NA, returnValue = NA)
path2ppi |
An object from the class |
species |
Either a number between 1 and the number of stored reference species or a
character string with the taxonomy id. If no value for |
returnValue |
Character value indicating whether to return a value. "proteins": a data frame
containing the proteins associated with the pathway of interest in the
corresponding reference species. "interactions": a data frame containing all
processed, relevant and non-redundant interactions. "irefindex": a data frame
containing all relevant interactions in the raw irefindex format. Is only
reasonable if |
See description for returnValue
Oliver Philipp [email protected]
addReference
, removeReference
,
showInteraction
data(ai) #Load test data set ppi <- Path2PPI("Autophagy induction", "Podospora anserina", "5145") ppi <- addReference(ppi, "Homo sapiens", "9606", human.ai.proteins, human.ai.irefindex, pa2human.ai.homologs) ppi <- addReference(ppi, "Saccharomyces cerevisiae (S288c)", "559292", yeast.ai.proteins, yeast.ai.irefindex, pa2yeast.ai.homologs) #Get general information about each stored reference species showReferences(ppi) #Get general information about reference species with the taxonomy id "9606" showReferences(ppi, species="9606") #Get all proteins associated with the pathway of interest #and previously given by the user proteins <- showReferences(ppi, species="9606", returnValue="proteins") #Get all processed and non-redundant interactions previously #determined to be relevant for the pathway of interest interactions <- showReferences(ppi, species="9606", returnValue="interactions") #Get all relevant interactions in the detailed irefindex format irefindex <- showReferences(ppi, species="9606", returnValue="irefindex")
data(ai) #Load test data set ppi <- Path2PPI("Autophagy induction", "Podospora anserina", "5145") ppi <- addReference(ppi, "Homo sapiens", "9606", human.ai.proteins, human.ai.irefindex, pa2human.ai.homologs) ppi <- addReference(ppi, "Saccharomyces cerevisiae (S288c)", "559292", yeast.ai.proteins, yeast.ai.irefindex, pa2yeast.ai.homologs) #Get general information about each stored reference species showReferences(ppi) #Get general information about reference species with the taxonomy id "9606" showReferences(ppi, species="9606") #Get all proteins associated with the pathway of interest #and previously given by the user proteins <- showReferences(ppi, species="9606", returnValue="proteins") #Get all processed and non-redundant interactions previously #determined to be relevant for the pathway of interest interactions <- showReferences(ppi, species="9606", returnValue="interactions") #Get all relevant interactions in the detailed irefindex format irefindex <- showReferences(ppi, species="9606", returnValue="irefindex")