2020-09-18 11:56:30 +00:00
|
|
|
# tocID <- "BIN-PPI-Analysis.R"
|
|
|
|
#
|
2017-09-12 20:09:20 +00:00
|
|
|
#
|
|
|
|
# Purpose: A Bioinformatics Course:
|
|
|
|
# R code accompanying the BIN-PPI-Analysis unit.
|
|
|
|
#
|
2020-10-04 04:11:13 +00:00
|
|
|
# Version: 1.3
|
2017-09-12 20:09:20 +00:00
|
|
|
#
|
2020-10-04 04:11:13 +00:00
|
|
|
# Date: 2017-08 - 2020-10
|
2017-09-12 20:09:20 +00:00
|
|
|
# Author: Boris Steipe (boris.steipe@utoronto.ca)
|
|
|
|
#
|
|
|
|
# Versions:
|
2020-10-04 04:11:13 +00:00
|
|
|
# 1.3 Bugfix: called the wrong function on ENSPsel in l. 220
|
2020-09-23 10:47:26 +00:00
|
|
|
# 1.2 2020 Updates; Rewrite for new STRINg V11;
|
|
|
|
# Deprecate save()/load() for saveRDS()/readRDS()
|
2019-01-08 07:11:25 +00:00
|
|
|
# 1.1 Change from require() to requireNamespace(),
|
|
|
|
# use <package>::<function>() idiom throughout,
|
|
|
|
# use Biocmanager:: not biocLite()
|
2017-11-10 04:44:25 +00:00
|
|
|
# 1.0 First live version
|
2017-09-12 20:09:20 +00:00
|
|
|
# 0.1 First code copied from 2016 material.
|
2017-11-10 04:44:25 +00:00
|
|
|
#
|
2017-09-12 20:09:20 +00:00
|
|
|
#
|
|
|
|
# TODO:
|
|
|
|
#
|
|
|
|
#
|
|
|
|
# == DO NOT SIMPLY source() THIS FILE! =======================================
|
2017-11-10 04:44:25 +00:00
|
|
|
#
|
2017-09-12 20:09:20 +00:00
|
|
|
# If there are portions you don't understand, use R's help system, Google for an
|
|
|
|
# answer, or ask your instructor. Don't continue if you don't understand what's
|
|
|
|
# going on. That's not how it works ...
|
2017-11-10 04:44:25 +00:00
|
|
|
#
|
2017-09-12 20:09:20 +00:00
|
|
|
# ==============================================================================
|
|
|
|
|
|
|
|
|
2017-11-10 04:44:25 +00:00
|
|
|
#TOC> ==========================================================================
|
2020-10-04 04:11:13 +00:00
|
|
|
#TOC>
|
2019-01-07 06:17:23 +00:00
|
|
|
#TOC> Section Title Line
|
|
|
|
#TOC> ---------------------------------------------------------------
|
2020-09-23 10:47:26 +00:00
|
|
|
#TOC> 1 Setup and data 49
|
|
|
|
#TOC> 2 Functional Edges in the Human Proteome 85
|
|
|
|
#TOC> 2.1 Cliques 128
|
|
|
|
#TOC> 2.2 Communities 169
|
|
|
|
#TOC> 2.3 Betweenness Centrality 183
|
|
|
|
#TOC> 3 biomaRt 230
|
|
|
|
#TOC> 4 Task for submission 301
|
2020-10-04 04:11:13 +00:00
|
|
|
#TOC>
|
2017-11-10 04:44:25 +00:00
|
|
|
#TOC> ==========================================================================
|
2017-09-12 20:09:20 +00:00
|
|
|
|
|
|
|
|
2017-11-10 04:44:25 +00:00
|
|
|
# = 1 Setup and data ======================================================
|
2017-09-12 20:09:20 +00:00
|
|
|
|
2017-11-10 04:44:25 +00:00
|
|
|
|
|
|
|
# Not surprisingly, the analysis of PPI networks needs iGraph:
|
|
|
|
|
2019-01-08 07:11:25 +00:00
|
|
|
if (! requireNamespace("igraph", quietly = TRUE)) {
|
2017-11-10 04:44:25 +00:00
|
|
|
install.packages("igraph")
|
|
|
|
}
|
|
|
|
# Package information:
|
|
|
|
# library(help = igraph) # basic information
|
|
|
|
# browseVignettes("igraph") # available vignettes
|
|
|
|
# data(package = "igraph") # available datasets
|
2017-09-12 20:09:20 +00:00
|
|
|
|
|
|
|
# In order for you to explore some real, biological networks, I give you a
|
2017-11-10 04:44:25 +00:00
|
|
|
# dataframe of functional relationships of human proteins that I have downloaded
|
|
|
|
# from the STRING database. The full table has 8.5 million records, here is a
|
|
|
|
# subset of records with combined confidence scores > 980
|
2017-09-12 20:09:20 +00:00
|
|
|
|
2020-09-23 10:47:26 +00:00
|
|
|
# The selected set of edges with a confidence of > 964 is a dataframe with about
|
|
|
|
# 50,000 edges and 8,400 unique proteins. Incidentaly, that's about the size of
|
2017-11-10 04:44:25 +00:00
|
|
|
# a fungal proteome. You can load the saved dataframe here (To read more about
|
2020-09-23 10:47:26 +00:00
|
|
|
# what the scores mean, see http://www.ncbi.nlm.nih.gov/pubmed/15608232 ).
|
2017-09-12 20:09:20 +00:00
|
|
|
|
2020-09-21 04:28:24 +00:00
|
|
|
STRINGedges <- readRDS("./data/STRINGedges.rds")
|
2017-09-12 20:09:20 +00:00
|
|
|
|
|
|
|
head(STRINGedges)
|
|
|
|
|
2017-11-10 04:44:25 +00:00
|
|
|
# Note that STRING has appended the tax-ID for Homo sapiens - 9606 - to the
|
|
|
|
# Ensemble transcript identifiers that start with ENSP. We'll remove them:
|
2017-09-12 20:09:20 +00:00
|
|
|
|
2020-09-23 10:47:26 +00:00
|
|
|
STRINGedges$a <- gsub("^9606\\.", "", STRINGedges$a)
|
|
|
|
STRINGedges$b <- gsub("^9606\\.", "", STRINGedges$b)
|
2017-11-10 04:44:25 +00:00
|
|
|
|
|
|
|
head(STRINGedges)
|
|
|
|
|
|
|
|
|
|
|
|
# = 2 Functional Edges in the Human Proteome ==============================
|
|
|
|
|
|
|
|
|
|
|
|
# There are many possibilities to explore interesting aspects of biological
|
|
|
|
# networks, we will keep with some very simple procedures here but you have
|
2020-09-23 10:47:26 +00:00
|
|
|
# to be aware that this is barely scratching the surface of possibilities.
|
2017-11-10 04:44:25 +00:00
|
|
|
# However, once the network exists in your computer, it is comparatively
|
2020-09-23 10:47:26 +00:00
|
|
|
# easy to find information online about the many, many options to analyze.
|
2017-11-10 04:44:25 +00:00
|
|
|
|
|
|
|
|
|
|
|
# Make a graph from this dataframe
|
2019-01-08 07:11:25 +00:00
|
|
|
?igraph::graph_from_data_frame
|
2017-09-12 20:09:20 +00:00
|
|
|
|
2019-01-08 07:11:25 +00:00
|
|
|
gSTR <- igraph::graph_from_data_frame(STRINGedges, directed = FALSE)
|
2017-09-12 20:09:20 +00:00
|
|
|
|
2020-09-23 10:47:26 +00:00
|
|
|
# CAUTION you DON'T want to plot a graph with 8,000 nodes and 50,000 edges -
|
2017-09-12 20:09:20 +00:00
|
|
|
# layout of such large graphs is possible, but requires specialized code. Google
|
|
|
|
# for <layout large graphs> if you are curious. Also, consider what one can
|
|
|
|
# really learn from plotting such a graph ...
|
|
|
|
|
|
|
|
# Of course simple computations on this graph are reasonably fast:
|
|
|
|
|
2019-01-08 07:11:25 +00:00
|
|
|
compSTR <- igraph::components(gSTR)
|
2017-09-12 20:09:20 +00:00
|
|
|
summary(compSTR) # our graph is fully connected!
|
|
|
|
|
2019-01-08 07:11:25 +00:00
|
|
|
hist(log(igraph::degree(gSTR)), col="#FEE0AF")
|
2017-09-12 20:09:20 +00:00
|
|
|
# this actually does look rather scale-free
|
|
|
|
|
2019-01-08 07:11:25 +00:00
|
|
|
(freqRank <- table(igraph::degree(gSTR)))
|
2017-09-12 20:09:20 +00:00
|
|
|
plot(log10(as.numeric(names(freqRank)) + 1),
|
|
|
|
log10(as.numeric(freqRank)), type = "b",
|
|
|
|
pch = 21, bg = "#FEE0AF",
|
|
|
|
xlab = "log(Rank)", ylab = "log(frequency)",
|
2020-09-23 10:47:26 +00:00
|
|
|
main = "8,400 nodes from the human functional interaction network")
|
2017-09-12 20:09:20 +00:00
|
|
|
|
|
|
|
# This looks very scale-free indeed.
|
2017-11-10 04:44:25 +00:00
|
|
|
|
|
|
|
(regressionLine <- lm(log10(as.numeric(freqRank)) ~
|
|
|
|
log10(as.numeric(names(freqRank)) + 1)))
|
|
|
|
abline(regressionLine, col = "firebrick")
|
|
|
|
|
2017-09-12 20:09:20 +00:00
|
|
|
# Now explore some more:
|
|
|
|
|
2017-11-10 04:44:25 +00:00
|
|
|
# == 2.1 Cliques ===========================================================
|
|
|
|
|
2017-09-12 20:09:20 +00:00
|
|
|
# Let's find the largest cliques. Remember: a clique is a fully connected
|
|
|
|
# subgraph, i.e. a subgraph in which every node is connected to every other.
|
|
|
|
# Biological complexes often appear as cliques in interaction graphs.
|
|
|
|
|
2019-01-08 07:11:25 +00:00
|
|
|
igraph::clique_num(gSTR)
|
2020-09-23 10:47:26 +00:00
|
|
|
# The largest clique has 81 members.
|
2017-09-12 20:09:20 +00:00
|
|
|
|
2019-01-08 07:11:25 +00:00
|
|
|
(C <- igraph::largest_cliques(gSTR)[[1]])
|
2017-09-12 20:09:20 +00:00
|
|
|
|
2020-09-23 10:47:26 +00:00
|
|
|
# Pick one of the proteins and find out what this fully connected cluster of 81
|
2017-11-10 04:44:25 +00:00
|
|
|
# proteins is (you can simply Google for any of the IDs). Is this expected?
|
|
|
|
|
|
|
|
# Plot this ...
|
2019-01-08 07:11:25 +00:00
|
|
|
R <- igraph::induced_subgraph(gSTR, C) # a graph from a selected set of vertices
|
2017-11-10 04:44:25 +00:00
|
|
|
|
|
|
|
# color the vertices along a color spectrum
|
2019-01-08 07:11:25 +00:00
|
|
|
vCol <- rainbow(igraph::gorder(R)) # "order" of a graph == number of nodes
|
2017-11-10 04:44:25 +00:00
|
|
|
|
|
|
|
# color the edges to have the same color as the originating node
|
|
|
|
eCol <- character()
|
|
|
|
for (i in seq_along(vCol)) {
|
2019-01-08 07:11:25 +00:00
|
|
|
eCol <- c(eCol, rep(vCol[i], igraph::gorder(R)))
|
2017-11-10 04:44:25 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
oPar <- par(mar= rep(0,4)) # Turn margins off
|
|
|
|
plot(R,
|
2019-01-08 07:11:25 +00:00
|
|
|
layout = igraph::layout_in_circle(R),
|
2017-11-10 04:44:25 +00:00
|
|
|
vertex.size = 3,
|
|
|
|
vertex.color = vCol,
|
|
|
|
edge.color = eCol,
|
|
|
|
edge.width = 0.1,
|
|
|
|
vertex.label = NA)
|
|
|
|
par(oPar)
|
2017-09-12 20:09:20 +00:00
|
|
|
|
2017-11-10 04:44:25 +00:00
|
|
|
# ... well: remember: a clique means every node is connected to every other
|
2020-09-23 10:47:26 +00:00
|
|
|
# node. We have 81 * 81 = 6,561 edges. This is what a matrix model of PPI
|
2017-11-10 04:44:25 +00:00
|
|
|
# networks looks like for large complexes.
|
2017-09-12 20:09:20 +00:00
|
|
|
|
|
|
|
|
2017-11-10 04:44:25 +00:00
|
|
|
# == 2.2 Communities =======================================================
|
|
|
|
|
2019-01-07 06:17:23 +00:00
|
|
|
set.seed(112358) # set RNG seed for repeatable randomness
|
2019-01-08 07:11:25 +00:00
|
|
|
gSTRclusters <- igraph::cluster_infomap(gSTR)
|
2019-01-07 06:17:23 +00:00
|
|
|
set.seed(NULL) # reset the RNG
|
|
|
|
|
2019-01-08 07:11:25 +00:00
|
|
|
igraph::modularity(gSTRclusters) # ... measures how separated the different
|
|
|
|
# membership types are from each other
|
|
|
|
tMem <- table(igraph::membership(gSTRclusters))
|
2020-09-23 10:47:26 +00:00
|
|
|
length(tMem) # About 700 communities identified
|
2019-01-08 07:11:25 +00:00
|
|
|
hist(tMem, breaks = 50, col = "skyblue") # most clusters are small ...
|
2020-09-23 10:47:26 +00:00
|
|
|
range(tMem) # ... but one has > 200 members
|
2017-11-10 04:44:25 +00:00
|
|
|
|
|
|
|
|
|
|
|
# == 2.3 Betweenness Centrality ============================================
|
2017-09-12 20:09:20 +00:00
|
|
|
|
|
|
|
# Let's find the nodes with the 10 - highest betweenness centralities.
|
|
|
|
#
|
2019-01-08 07:11:25 +00:00
|
|
|
BC <- igraph::centr_betw(gSTR)
|
2017-09-12 20:09:20 +00:00
|
|
|
|
|
|
|
# remember: BC$res contains the results
|
|
|
|
head(BC$res)
|
|
|
|
|
|
|
|
BC$res[1] # betweeness centrality of node 1 in the graph ...
|
|
|
|
# ... which one is node 1?
|
2019-01-08 07:11:25 +00:00
|
|
|
igraph::V(gSTR)[1]
|
2017-09-12 20:09:20 +00:00
|
|
|
|
2017-11-10 04:44:25 +00:00
|
|
|
# to get the ten-highest nodes, we simply label the elements of BC with their
|
2017-09-12 20:09:20 +00:00
|
|
|
# index ...
|
|
|
|
names(BC$res) <- as.character(1:length(BC$res))
|
|
|
|
|
|
|
|
# ... and then we sort:
|
|
|
|
sBC <- sort(BC$res, decreasing = TRUE)
|
|
|
|
head(sBC)
|
|
|
|
|
|
|
|
# This ordered vector means: node 3,862 has the highest betweeness centrality,
|
|
|
|
# node 1,720 has the second highest.. etc.
|
|
|
|
|
2017-11-10 04:44:25 +00:00
|
|
|
(BCsel <- as.numeric(names(sBC)[1:10]))
|
|
|
|
|
2017-09-12 20:09:20 +00:00
|
|
|
# We can use the first ten labels to subset the nodes in gSTR and fetch the
|
|
|
|
# IDs...
|
2019-01-08 07:11:25 +00:00
|
|
|
(ENSPsel <- names(igraph::V(gSTR)[BCsel]))
|
2017-09-12 20:09:20 +00:00
|
|
|
|
2020-09-23 10:47:26 +00:00
|
|
|
# Task:
|
|
|
|
# =====
|
|
|
|
# IMPORTANT, IF YOU INTEND TO SUBMIT YOUR ANALYSIS FOR CREDIT
|
2017-11-10 04:44:25 +00:00
|
|
|
# We are going to use these IDs to produce some output for a submitted task:
|
2020-09-23 10:47:26 +00:00
|
|
|
# therefore I need you to execute the following line, note the "seal" that this
|
2020-10-04 04:11:13 +00:00
|
|
|
# returns, and not change myENSPsel later:
|
2019-01-07 06:17:23 +00:00
|
|
|
|
2020-10-04 04:11:13 +00:00
|
|
|
myENSPsel <- selectENSP(ENSPsel)
|
2017-09-12 20:09:20 +00:00
|
|
|
|
|
|
|
# Next, to find what these proteins are...
|
|
|
|
|
|
|
|
# We could now Google for all of these IDs to learn more about them. But really,
|
|
|
|
# googling for IDs one after the other, that would be lame. Let's instead use
|
|
|
|
# the very, very useful biomaRt package to translate these Ensemble IDs into
|
|
|
|
# gene symbols.
|
|
|
|
|
2017-11-10 04:44:25 +00:00
|
|
|
|
|
|
|
# = 3 biomaRt =============================================================
|
|
|
|
|
2017-09-12 20:09:20 +00:00
|
|
|
|
|
|
|
# IDs are just labels, but for _bio_informatics we need to learn more about the
|
|
|
|
# biological function of the genes or proteins that we retrieve via graph data
|
|
|
|
# mining. biomaRt is the tool of choice. It's a package distributed by the
|
|
|
|
# bioconductor project. This here is not a biomaRt tutorial (that's for another
|
|
|
|
# day), simply a few lines of sample code to get you started on the specific use
|
|
|
|
# case of retrieving descriptions for ensembl protein IDs.
|
|
|
|
|
2019-01-08 07:11:25 +00:00
|
|
|
if (! requireNamespace("BiocManager", quietly = TRUE)) {
|
|
|
|
install.packages("BiocManager")
|
|
|
|
}
|
|
|
|
if (! requireNamespace("biomaRt", quietly = TRUE)) {
|
|
|
|
BiocManager::install("biomaRt")
|
2017-09-12 20:09:20 +00:00
|
|
|
}
|
2017-10-29 03:05:53 +00:00
|
|
|
# Package information:
|
|
|
|
# library(help = biomaRt) # basic information
|
|
|
|
# browseVignettes("biomaRt") # available vignettes
|
|
|
|
# data(package = "biomaRt") # available datasets
|
2017-09-12 20:09:20 +00:00
|
|
|
|
2020-09-23 10:47:26 +00:00
|
|
|
# define which dataset to use ... this takes a while for download
|
2019-01-08 07:11:25 +00:00
|
|
|
myMart <- biomaRt::useMart("ensembl", dataset="hsapiens_gene_ensembl")
|
2017-09-12 20:09:20 +00:00
|
|
|
|
|
|
|
# what filters are defined?
|
2020-09-23 10:47:26 +00:00
|
|
|
( filters <- biomaRt::listFilters(myMart) )
|
2017-11-10 04:44:25 +00:00
|
|
|
|
2017-09-12 20:09:20 +00:00
|
|
|
|
|
|
|
# and what attributes can we filter for?
|
2020-09-23 10:47:26 +00:00
|
|
|
( attributes <- biomaRt::listAttributes(myMart) )
|
2017-11-10 04:44:25 +00:00
|
|
|
|
2017-09-12 20:09:20 +00:00
|
|
|
|
|
|
|
# Soooo many options - let's look for the correct name of filters that are
|
|
|
|
# useful for ENSP IDs ...
|
|
|
|
filters[grep("ENSP", filters$description), ]
|
|
|
|
|
|
|
|
# ... and the correct attribute names for gene symbols and descriptions ...
|
2019-01-08 07:11:25 +00:00
|
|
|
attributes[grep("symbol", attributes$description, ignore.case = TRUE), ]
|
|
|
|
attributes[grep("description", attributes$description, ignore.case = TRUE), ]
|
2017-09-12 20:09:20 +00:00
|
|
|
|
|
|
|
|
|
|
|
# ... so we can put this together: here is a syntax example:
|
2019-01-08 07:11:25 +00:00
|
|
|
biomaRt::getBM(filters = "ensembl_peptide_id",
|
|
|
|
attributes = c("hgnc_symbol",
|
|
|
|
"wikigene_description",
|
|
|
|
"interpro_description",
|
|
|
|
"phenotype_description"),
|
|
|
|
values = "ENSP00000000442",
|
|
|
|
mart = myMart)
|
2017-09-12 20:09:20 +00:00
|
|
|
|
|
|
|
# A simple loop will now get us the information for our 10 most central genes
|
|
|
|
# from the human subset of STRING.
|
|
|
|
|
|
|
|
CPdefs <- list() # Since we don't know how many matches one of our queries
|
|
|
|
# will return, we'll put the result dataframes into a list.
|
|
|
|
|
2020-10-04 04:11:13 +00:00
|
|
|
for (ID in myENSPsel) {
|
2019-01-08 07:11:25 +00:00
|
|
|
CPdefs[[ID]] <- biomaRt::getBM(filters = "ensembl_peptide_id",
|
|
|
|
attributes = c("hgnc_symbol",
|
|
|
|
"wikigene_description",
|
|
|
|
"interpro_description",
|
|
|
|
"phenotype_description"),
|
|
|
|
values = ID,
|
|
|
|
mart = myMart)
|
2017-09-12 20:09:20 +00:00
|
|
|
}
|
|
|
|
|
2020-09-23 10:47:26 +00:00
|
|
|
|
2017-09-12 20:09:20 +00:00
|
|
|
# So what are the proteins with the ten highest betweenness centralities?
|
|
|
|
# ... are you surprised? (I am! Really.)
|
|
|
|
|
2017-11-10 04:44:25 +00:00
|
|
|
|
|
|
|
# = 4 Task for submission =================================================
|
|
|
|
|
|
|
|
# Write a loop that will go through your personalized list of Ensemble IDs and
|
2017-09-12 20:09:20 +00:00
|
|
|
# for each ID:
|
|
|
|
# -- print the ID,
|
2017-11-10 04:44:25 +00:00
|
|
|
# -- print the first row's HGNC symbol,
|
2017-09-12 20:09:20 +00:00
|
|
|
# -- print the first row's wikigene description.
|
2017-11-10 04:44:25 +00:00
|
|
|
# -- print the first row's phenotype.
|
2017-09-12 20:09:20 +00:00
|
|
|
#
|
2020-09-23 10:47:26 +00:00
|
|
|
# Write your thoughts about this group of genes.
|
|
|
|
#
|
2017-09-12 20:09:20 +00:00
|
|
|
# (Hint, you can structure your loop in the same way as the loop that
|
|
|
|
# created CPdefs. )
|
|
|
|
|
2020-09-23 10:47:26 +00:00
|
|
|
# Submit the "seal" for your ENSP vector, the ENSP vector itself, the R code
|
|
|
|
# for this loop and its output into your report if you are submitting
|
|
|
|
# anything for credit for this unit. Please read the requirements carefully.
|
2017-09-12 20:09:20 +00:00
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
# [END]
|