Gene Ontology analysis with Python and Bioconductor

Brad Chapman bio photo By Brad Chapman Comment

Last post we discussed extracting differentially expressed genes from a set of short read data. The output of this analysis is a large list of genes, and the next challenge is to try and extract interesting patterns from the results. Are genes in particular families over-represented? Are there a large number of genes with a particular function? These patterns can help direct your future experiments and provide an understandable grouping for biologists examining the data.

The Gene Ontology (GO) project provides a standardized set of terms describing the molecular function of genes. We will use the topGO package from the Bioconductor project to identify over-represented GO terms from a set of differentially expressed genes. Python will be used to prepare the data, utilizing rpy2 to call R for the statistical analysis.

The first step is to parse input files describing the differentially expressed genes and the mapping of gene names to GO terms. For the example we will use a simple CSV file from our previous analysis and an equally simple file describing the gene to GO mapping. For real life work, you will need to get the GO mapping file for your organism; for instance, here is the Arabidopsis GO mapping from TAIR.

The input file of differentially expressed genes is parsed into a dictionary where the keys are gene names and the values are associated p-values:

[sourcecode language="python"]
import csv

def parse_input_csv(in_handle):
reader = csv.reader(in_handle)
reader.next() # header
all_genes = dict()
for (gene_name, _, _, pval) in reader:
all_genes[gene_name] = float(pval)
return all_genes
[/sourcecode]

The GO mapping input file is parsed into two dictionaries: one mapping genes to GO terms for the GO analysis, and the second mapping GO terms to genes for displaying the results of the analysis:

[sourcecode language="python"]
import collections

def parse_go_map_file(in_handle, genes_w_pvals):
gene_list = genes_w_pvals.keys()
gene_to_go = collections.defaultdict(list)
go_to_gene = collections.defaultdict(list)
for line in in_handle:
parts = line.split("\t")
gene_id = parts[0]
go_id = parts[1].strip()
if gene_id in gene_list:
gene_to_go[gene_id].append(go_id)
go_to_gene[go_id].append(gene_id)
return dict(gene_to_go), dict(go_to_gene)
[/sourcecode]

Next we set up our associated R environment for GO analysis. This involves loading the topGO library and creating a function to classify genes as either differentially expressed, or not. We do this based on the p-value, classifying genes with a p-value of 0.01 or less as differentially expressed:

[sourcecode language="python"]
import rpy2.robjects as robjects

gene_pval = 1e-2
robjects.r('''
library(topGO)
''')
robjects.r('''
topDiffGenes = function(allScore) {
return (allScore < %s)
}
''' % gene_pval)
[/sourcecode]

Now it's time to run the analysis. The paramters are defined as a dictionary which initializes a topGOdata object:

  • ontology -- the type of GO ontology to analyze; molecular function (MF) is examined here.
  • annot -- how the GO terms are annotated; we specify they will be identified by a mapping of genes to GO ids, passed as the gene2GO argument.
  • geneSelectionFun -- the function we defined above to determine if genes are differentially expressed.
  • allGenes -- an R named vector where the names are genes and the values are p-values; we use a utility function to convert a python dictionary to the named vector.
  • gene2GO -- an R named vector mapping gene names to associated GO terms.

The initialized object is run using Fisher classic analysis. Several different analysis methods and test statistics are available in the package. The resulting scores are extracted and a summary table of the results is generated:

[sourcecode language="python"]
def _dict_to_namedvector(init_dict):
"""Call R to create a named vector from an input dictionary.
"""
return robjects.r.c(**init_dict)

go_term_type = "MF"
topgo_method = "classic" # choice of classic, elim, weight
params = {"ontology" : go_term_type,
"annot" : robjects.r["annFUN.gene2GO"],
"geneSelectionFun" : robjects.r["topDiffGenes"],
"allGenes" : _dict_to_namedvector(gene_vals),
"gene2GO" : _dict_to_namedvector(gene_to_go)
}
go_data = robjects.r.new("topGOdata", **params)
results = robjects.r.runTest(go_data, algorithm=topgo_method,
statistic="fisher")
scores = robjects.r.score(results)
score_names = scores.getnames()
num_summarize = min(100, len(score_names))
results_table = robjects.r.GenTable(go_data, elimFisher=results,
orderBy="elimFisher", topNodes=num_summarize)
[/sourcecode]

The scores and results table are used to generate a list of over-represented GO terms with their associated p-values and descriptions:

[sourcecode language="python"]
# extract term names from the topGO summary dataframe
GO_ID_INDEX = 0
TERM_INDEX = 1
ids_to_terms = dict()
for index, go_id in enumerate(results_table[GO_ID_INDEX]):
ids_to_terms[go_id] = results_table[TERM_INDEX][index]
go_terms = []
# convert the scores and results information info terms to return
for index, item in enumerate(scores):
if item < go_pval:
go_id = score_names[index]
go_terms.append((item, go_id, ids_to_terms.get(go_id, "")))
go_terms.sort()
[/sourcecode]

Finally, we print out the resulting overexpressed GO terms along with associated genes:

[sourcecode language="python"]
def print_go_info(go_terms, go_term_type, go_to_gene):
for final_pval, go_id, go_term in go_terms:
genes = []
for check_go in [go_id] + get_go_children(go_id, go_term_type):
genes.extend(go_to_gene.get(check_go, []))
genes = sorted(list(set(genes)))
print "-> %s (%s) : %0.4f" % (go_id, go_term, final_pval)
for g in genes:
print g
[/sourcecode]

One tricky part of the post-analysis processing is finding all of the genes associated with an identified GO id. The GO terminology is hierarchical, so a resulting GO identifier like GO:0042578 in the example results file can represent both that term and more specific terms. A nice way to illustrate this for yourself is to look at the AmiGO browser display for the term.

The result is that genes associated with an over-represented term will come from both the identified term and more specific terms. Bioconductor has an interface to the GO database, so we can extract all of the relevant terms by moving down the GO tree collecting the more specific children, given the parent term:

[sourcecode language="python"]
def get_go_children(go_term, go_term_type):
robjects.r('''
library(GO.db)
''')
child_map = robjects.r["GO%sCHILDREN" % (go_term_type)]
children = []
to_check = [go_term]
while len(to_check) > 0:
new_children = []
for check_term in to_check:
new_children.extend(list(robjects.r.get(check_term, child_map)))
new_children = list(set([c for c in new_children if c]))
children.extend(new_children)
to_check = new_children
children = list(set(children))
return children
[/sourcecode]

The full script pulls all of this code together into a working example. You will need to update the parsers to match your input and GO mapping file, but otherwise the code can plug directly into full genome analyses of differential expression data.

comments powered by Disqus