Bioconductor has extensive facilities for mapping between microarray probe, gene, pathway, gene ontology, homology and other annotations.
Bioconductor has built-in representations of GO, KEGG, vendor, and other annotations, and can easily access NCBI, Biomart, UCSC, and other sources.
The following illustrates a typical R / Bioconductor session for a ChipDb package. It continues the differential expression workflow, taking a 'top table' of differentially expressed probesets and discovering the genes probed, and the Gene Ontology pathways to which they belong.
First lets consider some typical probe Ids. If you have done a microarray analysis before you have probably already run into IDs like this. They are typically manufacturer assigned and normally only relevant to a small number of chips. Below I am just going to demonstrate on 6 probe Ids from the u133 2.0 affymetrix platform.
## Affymetrix U133 2.0 array IDs of interest; these might be
## obtained from
##
## tbl <- topTable(efit, coef=2)
## ids <- tbl[["ID"]]
##
## as part of a more extensive workflow.
ids <- c("39730_at", "1635_at", "1674_at", "40504_at", "40202_at")
Load libraries as sources of annotation
library("hgu95av2.db")
To list the kinds of things that can be retrieved, use the columns method.
columns(hgu95av2.db)
## [1] "PROBEID" "ENTREZID" "PFAM" "IPI"
## [5] "PROSITE" "ACCNUM" "ALIAS" "CHR"
## [9] "CHRLOC" "CHRLOCEND" "ENZYME" "MAP"
## [13] "PATH" "PMID" "REFSEQ" "SYMBOL"
## [17] "UNIGENE" "ENSEMBL" "ENSEMBLPROT" "ENSEMBLTRANS"
## [21] "GENENAME" "UNIPROT" "GO" "EVIDENCE"
## [25] "ONTOLOGY" "GOALL" "EVIDENCEALL" "ONTOLOGYALL"
## [29] "OMIM" "UCSCKG"
To list the kinds of things that can be used as keys we can use the keytypes method
keytypes(hgu95av2.db)
## [1] "ENTREZID" "PFAM" "IPI" "PROSITE"
## [5] "ACCNUM" "ALIAS" "CHR" "CHRLOC"
## [9] "CHRLOCEND" "ENZYME" "MAP" "PATH"
## [13] "PMID" "REFSEQ" "SYMBOL" "UNIGENE"
## [17] "ENSEMBL" "ENSEMBLPROT" "ENSEMBLTRANS" "GENENAME"
## [21] "UNIPROT" "GO" "EVIDENCE" "ONTOLOGY"
## [25] "GOALL" "EVIDENCEALL" "ONTOLOGYALL" "PROBEID"
## [29] "OMIM" "UCSCKG"
To extract viable keys of a particular kind, use the keys method.
head(keys(hgu95av2.db, keytype="ENTREZID"))
## [1] "10" "100" "1000" "10000" "100008586" "10001"
The select method allows you to map probe ids to ENTREZ gene ids…
select(hgu95av2.db, keys=ids, columns="ENTREZID", keytype="PROBEID")
## PROBEID ENTREZID
## 1 39730_at 25
## 2 1635_at 25
## 3 1674_at 7525
## 4 40504_at 5445
## 5 40202_at 687
Or to gene names and gene symbols etc.
select(hgu95av2.db, keys=ids, columns=c("ENTREZID","GENENAME", "SYMBOL"), keytype="PROBEID")
## PROBEID ENTREZID GENENAME
## 1 39730_at 25 c-abl oncogene 1, non-receptor tyrosine kinase
## 2 1635_at 25 c-abl oncogene 1, non-receptor tyrosine kinase
## 3 1674_at 7525 v-yes-1 Yamaguchi sarcoma viral oncogene homolog 1
## 4 40504_at 5445 paraoxonase 2
## 5 40202_at 687 Kruppel-like factor 9
## SYMBOL
## 1 ABL1
## 2 ABL1
## 3 YES1
## 4 PON2
## 5 KLF9
We can also find and extract the GO ids associated with the first id (there are quite a few of these)
res <- select(hgu95av2.db, keys=ids[1], columns="GO", keytype="PROBEID")
## Warning: 'select' resulted in 1:many mapping between keys and return rows
head(res)
## PROBEID GO EVIDENCE ONTOLOGY
## 1 39730_at GO:0000287 IDA MF
## 2 39730_at GO:0003677 NAS MF
## 3 39730_at GO:0003785 TAS MF
## 4 39730_at GO:0004515 TAS MF
## 5 39730_at GO:0004713 IDA MF
## 6 39730_at GO:0004715 IDA MF
And then once we have done that, we can also use the GO.db package to find the Terms associated with those GOIDs. How will this work? Well the GO.db package will load a GO.db object that can be used in a manner similar to what we just saw with our ChipDb object hgu95av2.db. So we can use the same four methods that we just learned about (columns, keytypes, keys and select), to extract whatever data we need from this object.
library("GO.db")
##
head(res$GO) ## shows what we are using as keys
## [1] "GO:0000287" "GO:0003677" "GO:0003785" "GO:0004515" "GO:0004713"
## [6] "GO:0004715"
head(select(GO.db, keys=res$GO, columns="TERM", keytype="GOID"))
## GOID TERM
## 1 GO:0000287 magnesium ion binding
## 2 GO:0003677 DNA binding
## 3 GO:0003785 actin monomer binding
## 4 GO:0004515 nicotinate-nucleotide adenylyltransferase activity
## 5 GO:0004713 protein tyrosine kinase activity
## 6 GO:0004715 non-membrane spanning protein tyrosine kinase activity
[ Back to top ]
The organism wide gene centered packages (OrgDb packages) all contain gene centered data for an organism. These packages are accessed in much the same way as the platform based (ChipDb) package and the Go.db package the we just discussed. But because they are general, they don't contain information like probe IDs that would relate to any specific microarray platform.
But the more important thing to understand is that the same methods apply. So for example you can look up information in this way:
library(org.Hs.eg.db)
keys <- head(keys(org.Hs.eg.db, keytype="ENTREZID"), n=2)
columns <- c("PFAM","GO", "SYMBOL")
select(org.Hs.eg.db, keys, columns, keytype="ENTREZID")
## Warning: 'select' resulted in 1:many mapping between keys and return rows
## ENTREZID PFAM GO EVIDENCE ONTOLOGY SYMBOL
## 1 1 <NA> GO:0003674 ND MF A1BG
## 2 1 <NA> GO:0005576 IDA CC A1BG
## 3 1 <NA> GO:0008150 ND BP A1BG
## 4 1 <NA> GO:0070062 IDA CC A1BG
## 5 1 <NA> GO:0072562 IDA CC A1BG
## 6 10 PF00797 GO:0004060 IEA MF NAT2
## 7 10 PF00797 GO:0005829 TAS CC NAT2
## 8 10 PF00797 GO:0006805 TAS BP NAT2
## 9 10 PF00797 GO:0044281 TAS BP NAT2
[ Back to top ]
Exercise 1: Look at the help page for the different columns and keytypes values with: help(“SYMBOL”). Now use this information and what we just described to look up the entrez gene, probe id and chromosome for the gene symbol “MSX2”.
Exercise 2: Examine the gene symbols for both the hgu95av2.db and the org.Hs.eg.db packages. Which one has more gene symbols? Which one has more gene symbols that can be mapped to an entrez gene ID? Which object seems to contain more information?
Exercise 3: In the previous exercise we had to use gene symbols as keys. But in the past this kind of behavior has sometimes been inadvisable because some gene symbols are used as the official symbol for more than one gene. To learn if this is still happening take advantage of the fact that entrez gene ids are uniquely assigned, and extract all of the gene symbols and their associated entrez gene ids from the org.Hs.eg.db package. Then check the symbols for redundancy.
[ Back to top ]
The genome centered TranscriptDb packages support the same interface as that ChipDb and the OrgDb packages.
library(TxDb.Hsapiens.UCSC.hg19.knownGene)
txdb <- TxDb.Hsapiens.UCSC.hg19.knownGene ## done for convenience
keys <- head(keys(txdb, keytype="GENEID"), n=2)
columns <- c("TXNAME", "TXSTART","TXSTRAND")
select(txdb, keys, columns, keytype="GENEID")
## Warning: 'select' resulted in 1:many mapping between keys and return rows
## GENEID TXNAME TXSTRAND TXSTART
## 1 1 uc002qsd.4 - 58858172
## 2 1 uc002qsf.2 - 58859832
## 3 10 uc003wyw.1 + 18248755
But in addition to supporting the standard set of methods (select, keytypes, keys and columns). The TranscriptDb objects also support methods to retrieve the annotations as ranges. These accessors break down into two basic categories. The most basic will return annotations as GRanges objects. Some examples of these are: transcripts(), exons() and cds().
This for example will return all the transcripts as ranges:
transcripts(txdb)
## GRanges with 82960 ranges and 2 metadata columns:
## seqnames ranges strand | tx_id tx_name
## <Rle> <IRanges> <Rle> | <integer> <character>
## [1] chr1 [ 11874, 14409] + | 1 uc001aaa.3
## [2] chr1 [ 11874, 14409] + | 2 uc010nxq.1
## [3] chr1 [ 11874, 14409] + | 3 uc010nxr.1
## [4] chr1 [ 69091, 70008] + | 4 uc001aal.1
## [5] chr1 [321084, 321115] + | 5 uc001aaq.2
## ... ... ... ... ... ... ...
## [82956] chrY [27605645, 27605678] - | 78803 uc004fwx.1
## [82957] chrY [27606394, 27606421] - | 78804 uc022cpc.1
## [82958] chrY [27607404, 27607432] - | 78805 uc004fwz.3
## [82959] chrY [27635919, 27635954] - | 78806 uc022cpd.1
## [82960] chrY [59358329, 59360854] - | 78807 uc011ncc.1
## ---
## seqlengths:
## chr1 chr2 ... chrUn_gl000249
## 249250621 243199373 ... 38502
And this will return all the exons as ranges:
exons(txdb)
## GRanges with 289969 ranges and 1 metadata column:
## seqnames ranges strand | exon_id
## <Rle> <IRanges> <Rle> | <integer>
## [1] chr1 [11874, 12227] + | 1
## [2] chr1 [12595, 12721] + | 2
## [3] chr1 [12613, 12721] + | 3
## [4] chr1 [12646, 12697] + | 4
## [5] chr1 [13221, 14409] + | 5
## ... ... ... ... ... ...
## [289965] chrY [27607404, 27607432] - | 277746
## [289966] chrY [27635919, 27635954] - | 277747
## [289967] chrY [59358329, 59359508] - | 277748
## [289968] chrY [59360007, 59360115] - | 277749
## [289969] chrY [59360501, 59360854] - | 277750
## ---
## seqlengths:
## chr1 chr2 ... chrUn_gl000249
## 249250621 243199373 ... 38502
But these operations will also support the extraction of extra metadata. All extra data will be inserted into the metadata slot of the returned GRanges object. So for example you could spice up your call to transcripts by using the columns argument like this.
transcripts(txdb, columns = c("tx_id","tx_name","gene_id"))
## GRanges with 82960 ranges and 3 metadata columns:
## seqnames ranges strand | tx_id tx_name
## <Rle> <IRanges> <Rle> | <integer> <character>
## [1] chr1 [ 11874, 14409] + | 1 uc001aaa.3
## [2] chr1 [ 11874, 14409] + | 2 uc010nxq.1
## [3] chr1 [ 11874, 14409] + | 3 uc010nxr.1
## [4] chr1 [ 69091, 70008] + | 4 uc001aal.1
## [5] chr1 [321084, 321115] + | 5 uc001aaq.2
## ... ... ... ... ... ... ...
## [82956] chrY [27605645, 27605678] - | 78803 uc004fwx.1
## [82957] chrY [27606394, 27606421] - | 78804 uc022cpc.1
## [82958] chrY [27607404, 27607432] - | 78805 uc004fwz.3
## [82959] chrY [27635919, 27635954] - | 78806 uc022cpd.1
## [82960] chrY [59358329, 59360854] - | 78807 uc011ncc.1
## gene_id
## <CharacterList>
## [1] 100287102
## [2] 100287102
## [3] 100287102
## [4] 79501
## [5]
## ... ...
## [82956]
## [82957]
## [82958]
## [82959]
## [82960]
## ---
## seqlengths:
## chr1 chr2 ... chrUn_gl000249
## 249250621 243199373 ... 38502
The 2nd kind of range accessor supported by TranscriptDb objects are the ones that return GRangesList objects. Some examples of these are: transcriptsBy(), exonsBy() or cdsBy(). These accessors just allow you to return a GRangesList object that contains the desired ranges by split up by some important feature type that is specified using the “by” argument. A typical case is to extract all the transcript ranges known for all the genes. You can do that like this:
transcriptsBy(txdb, by="gene")
## GRangesList of length 23459:
## $1
## GRanges with 2 ranges and 2 metadata columns:
## seqnames ranges strand | tx_id tx_name
## <Rle> <IRanges> <Rle> | <integer> <character>
## [1] chr19 [58858172, 58864865] - | 70455 uc002qsd.4
## [2] chr19 [58859832, 58874214] - | 70456 uc002qsf.2
##
## $10
## GRanges with 1 range and 2 metadata columns:
## seqnames ranges strand | tx_id tx_name
## [1] chr8 [18248755, 18258723] + | 31944 uc003wyw.1
##
## $100
## GRanges with 1 range and 2 metadata columns:
## seqnames ranges strand | tx_id tx_name
## [1] chr20 [43248163, 43280376] - | 72132 uc002xmj.3
##
## ...
## <23456 more elements>
## ---
## seqlengths:
## chr1 chr2 ... chrUn_gl000249
## 249250621 243199373 ... 38502
[ Back to top ]
Exercise 4: Use the accessors for the TxDb.Hsapiens.UCSC.hg19.knownGene package to retrieve the gene id, transcript name and transcript chromosome for all the transcripts. Do this using both the select() method and also using the transcripts() method. What is the difference in the output?
Exercise 5: Load the TxDb.Athaliana.BioMart.plantsmart16 package. This package is not from UCSC and it is based on plantsmart. Now use select or one of the range based accessors to look at the gene ids from this TranscriptDb object. How tdo they compare to what you saw in the TxDb.Hsapiens.UCSC.hg19.knownGene package?
[ Back to top ]
What if you wanted to combine all the good stuff from the GO.db package with what you find in the appropriate TranscriptDb and OrgDb packages for an organism? Then you would want to use an OrganismDb package. An example of an OrganismDb package is the Homo.sapiens package. Like the OrgDb, ChipDb and TranscriptDb packages, it supports the use of select, keytypes, keys and columns.
library(Homo.sapiens)
keys <- head(keys(Homo.sapiens, keytype="ENTREZID"), n=2)
columns <- c("SYMBOL","TXNAME")
select(Homo.sapiens, keys, columns, keytype="ENTREZID")
## Warning: 'select' resulted in 1:many mapping between keys and return rows
## Warning: 'select' resulted in 1:many mapping between keys and return rows
## ENTREZID SYMBOL TXNAME
## 1 1 A1BG uc002qsd.4
## 2 1 A1BG uc002qsf.2
## 3 10 NAT2 uc003wyw.1
When an OrganismDb package knows about a relevant TranscriptDb package, it can also support the ranged accessors introduced with the TranscriptDb objects.
transcripts(Homo.sapiens, columns=c("TXNAME","SYMBOL"))
## GRanges with 82960 ranges and 2 metadata columns:
## seqnames ranges strand | TXNAME
## <Rle> <IRanges> <Rle> | <CharacterList>
## [1] chr1 [ 11874, 14409] + | uc001aaa.3
## [2] chr1 [ 11874, 14409] + | uc010nxq.1
## [3] chr1 [ 11874, 14409] + | uc010nxr.1
## [4] chr1 [ 69091, 70008] + | uc001aal.1
## [5] chr1 [321084, 321115] + | uc001aaq.2
## ... ... ... ... ... ...
## [82956] chrY [27605645, 27605678] - | uc004fwx.1
## [82957] chrY [27606394, 27606421] - | uc022cpc.1
## [82958] chrY [27607404, 27607432] - | uc004fwz.3
## [82959] chrY [27635919, 27635954] - | uc022cpd.1
## [82960] chrY [59358329, 59360854] - | uc011ncc.1
## SYMBOL
## <CharacterList>
## [1] DDX11L1
## [2] DDX11L1
## [3] DDX11L1
## [4] OR4F5
## [5] NA
## ... ...
## [82956] NA
## [82957] NA
## [82958] NA
## [82959] NA
## [82960] NA
## ---
## seqlengths:
## chr1 chr2 ... chrUn_gl000249
## 249250621 243199373 ... 38502
You might be surprised to learn that an OrganismDb package does not itself contain very much information. Instead, it “knows where to find it”, but referencing other packages that themselves implement a select interface. So to create an OrganismDb package, you really only need to specify where the information needs to come from. Configuring an OrganismDb object is therefore pretty simple. You simply create a special list object that describes which IDs from each package are the same kind of IDs in other packages to be included, along with the relevant package names. So in the following example, the “GOID” values from the GO.db package act as foreign keys for the “GO” values in the org.Hs.eg.db package and so on.
gd <- list(join1 = c(GO.db="GOID", org.Hs.eg.db="GO"),
join2 = c(org.Hs.eg.db="ENTREZID",
TxDb.Hsapiens.UCSC.hg19.knownGene="GENEID"))
makeOrganismPackage(pkgname = "Homo.sapiens",
graphData = gd,
organism = "Homo sapiens",
version = "1.0.0",
maintainer = "Package Maintainer<maintainer@somewhere.org>",
author = "Some Body",
destDir = ".",
license = "Artistic-2.0")
In this way, you can create a custom OrganismDb package for any organism of interest, providing that you have also have access to the supporting packages. There is a vignette that covers this topic in more detail here.
Exercise 6: Use the Homo.sapiens object to look up the gene symbol, transcript start and chromosome using select(). Then do the same thing using transcripts. You might expect that this call to transcripts will look the same as it did for the TranscriptDb object, but (temporarily) it will not.
Exercise 7: Look at the results from call the columns method on the Homo.sapiens object and compare that to what happens when you call columns on the org.Hs.eg.db object and then look at a call to columns on the TxDb.Hsapiens.UCSC.hg19.knownGene object. What is the difference between TXSTART and CHRLOC? Which one do you think you should use for transcripts or other genomic information?
[ Back to top ]
Lets look more closely at the keys method. We have already talked about how you can use it to do this:
library(Homo.sapiens)
keys <- head(keys(Homo.sapiens, keytype="ENTREZID"), n=2)
And then you can use it with select to look up other kinds of information. But what if you only know partial information about the keys you are looking up? In Bioconductor 2.13 and higher there are extra arguments for the keys method that you can make use of to find keys that match certain criteria. The most useful is probably the pattern argument. The pattern argument allows you to find out which keys match a certain pattern. So for example, you can look up entrez gene IDs that start with a “2” like this:
head(keys(Homo.sapiens, keytype="ENTREZID", pattern="^2"), n=6)
## [1] "2" "20" "2000" "200008" "200010" "200014"
Or you could look up gene symbols that start with “MS”:
head(keys(Homo.sapiens, keytype="SYMBOL", pattern="^MS"), n=6)
## [1] "MS4A1" "MS4A3" "MS4A2" "MSTN" "MSH6" "MS"
If your string matching is too specific, you could also try to use fuzzy matching by setting the fuzzy argument to TRUE:
head(keys(Homo.sapiens, keytype="SYMBOL", pattern="^MS", fuzzy=TRUE), n=6)
## [1] "MS4A1" "MS4A3" "MS4A2" "MSTN" "MSH6" "LIMS1"
And if you want to match one one key and actually return another, then you can use the column argument to indicate which key you want to search for pattern on while using the keytype to indicate which kind of key you want returned. So you could (for example) get back ensembl IDs where the symbol starts with “MS”.
keys <- head(keys(Homo.sapiens, keytype="ENSEMBL", pattern="^MS", column="SYMBOL"), n=6)
keys
## [1] "ENSG00000156738" "ENSG00000149516" "ENSG00000138379" "ENSG00000116062"
## [5] "ENSG00000095002" "ENSG00000113318"
select(Homo.sapiens, keys, "SYMBOL", keytype="ENSEMBL")
## ENSEMBL SYMBOL
## 1 ENSG00000156738 MS4A1
## 2 ENSG00000149516 MS4A3
## 3 ENSG00000138379 MSTN
## 4 ENSG00000116062 MSH6
## 5 ENSG00000095002 MSH2
## 6 ENSG00000113318 MSH3
Exercise 8: Use the Homo.sapiens object with the keys method to look up the entrez gene IDs for all gene symbols that contain the letter “X”.
[ Back to top ]
So far we have been discussing annotations that are fairly well established and that represent consensus findings from the scientific community. These kinds of annotations are usually curated at large governmental institutions like NCBI or ensembl and for the most part everyone basically agrees about what they mean and how to use them.
But sometimes the annotations that you need are not as well established. Sometimes (for example) we just need to compare our results to the data from a recent large study such as the encode project. The AnnotationHub package is designed to be useful for getting access to data like this. AnnotationHub allows you to get access to data from a range of different data reposotories, with the caveat that the data objects in AnnotationHub have all been pre-processed into appropriate R objects for you.
To make use of AnnotationHub, you need to load the package and then create an AnnotationHub object. Notice that unlike the other packages, with AnnotationHub, you have to create an AnnotationHub object when you 1st start up your R session.
library(AnnotationHub)
ah = AnnotationHub()
Once you have done this, you can “find” any of the available resources just by tab completing along a path like this.
res <- ah$goldenpath.hg19.encodeDCC.wgEncodeUwTfbs.wgEncodeUwTfbsMcf7CtcfStdPkRep1.narrowPeak_0.0.1.RData
res
## GRanges with 82163 ranges and 6 metadata columns:
## seqnames ranges strand | name score
## <Rle> <IRanges> <Rle> | <character> <integer>
## [1] chr1 [237640, 237790] * | . 0
## [2] chr1 [544660, 544810] * | . 0
## [3] chr1 [567480, 567630] * | . 0
## [4] chr1 [569820, 569970] * | . 0
## [5] chr1 [714200, 714350] * | . 0
## ... ... ... ... ... ... ...
## [82159] chrX [154764540, 154764690] * | . 0
## [82160] chrX [154807400, 154807550] * | . 0
## [82161] chrX [154881060, 154881210] * | . 0
## [82162] chrX [154892100, 154892250] * | . 0
## [82163] chrX [154916040, 154916190] * | . 0
## signalValue pValue qValue peak
## <numeric> <numeric> <numeric> <integer>
## [1] 30 26.892 -1 -1
## [2] 6 8.164 -1 -1
## [3] 100 56.718 -1 -1
## [4] 85 49.654 -1 -1
## [5] 17 13.184 -1 -1
## ... ... ... ... ...
## [82159] 26 25.29 -1 -1
## [82160] 22 27.65 -1 -1
## [82161] 17 16.42 -1 -1
## [82162] 72 101.61 -1 -1
## [82163] 32 32.52 -1 -1
## ---
## seqlengths:
## chr1 chr10 chr11 ... chr8 chr9 chrX
## 249250621 135534747 135006516 ... 146364022 141213431 155270560
In the above example, AnnotationHub will retrieve, download and cache locally the file that you tab-completed to, and then store the results in “res”.
Now you can see how many ways there are to currently complete that path, by checking the length of the AnnotationHub object:
length(ah)
## [1] 10780
The AnnotationHub is still a pretty new resource, and we already hav a LOT of things in there! How can we narrow this down? Right now we can use filters. By default, there are no filters applied, so calling filters() on our AnnotationHub is just an empty list.
filters(ah)
## list()
What things can be used as filters? We can use the columns() method to find out.
columns(ah)
## [1] "BiocVersion" "DataProvider" "Title"
## [4] "SourceFile" "Species" "SourceUrl"
## [7] "SourceVersion" "TaxonomyId" "Genome"
## [10] "Description" "Tags" "RDataClass"
## [13] "RDataPath" "Coordinate_1_based" "Maintainer"
## [16] "RDataVersion" "RDataDateAdded" "Recipe"
What values can be used with these filters? Here, the keys method will give us an answer.
head(keys(ah, keytype="Species"))
## [1] "9606" "Acromyrmex echinatior"
## [3] "Acyrthosiphon pisum" "Aedes aegypti"
## [5] "Agaricus bisporus" "Ailuropoda melanoleuca"
So now we know what we need to apply a filter to our AnnotationHub. The following filter will limit our AnnotationHub to just those entries that correspond to cattle (Bos taurus).
filters(ah) <- list(Species="Bos taurus")
length(ah)
## [1] 145
We can also view and filter our AnnotationHub object interactively by simply calling the display function on it
d <- display(ah)
We can then filter the AnnotationHub object for “Homo sapiens” by either using the Global search field on the top right corner of the page or the in-column search field for “Species”.
By default 1000 entries are displayed per page, we can change this using the filter on the top of the page or navigate through different pages using the page scrolling feature at the bottom of the page.
We can also select the rows of interest to us and send them back to the R session using 'Send Rows' button ; this sets a filter internally which filters the AnnotationHub object.
[ Back to top ]
Exercise 9: Set the AnnotationHub filter to NULL to clear it out, and then set ip up so that it is extracting data that originated with the UCSC data provider and that is also limited to Homo sapiens and the hg19 genome.
Exercise 10 Now that you have basically narrowed things down to the hg19 annotations from UCSC genome browser, lets get one of these annotations. Now tab complete your way to the oreganno track and save it into a local variable.
[ Back to top ]
Another valuable resource is the biomaRt package. The biomaRt package exposes a huge family of online annotation resources called marts. Here is a brief run down of how to use it. For the first step, load the package and decide which “mart” you want to use, then use the useMart() method to create a mart object
library("biomaRt")
head(listMarts())
## biomart version
## 1 ensembl ENSEMBL GENES 75 (SANGER UK)
## 2 snp ENSEMBL VARIATION 75 (SANGER UK)
## 3 functional_genomics ENSEMBL REGULATION 75 (SANGER UK)
## 4 vega VEGA 53 (SANGER UK)
## 5 fungi_mart_22 ENSEMBL FUNGI 22 (EBI UK)
## 6 fungi_variations_22 ENSEMBL FUNGI VARIATION 22 (EBI UK)
ensembl <- useMart("ensembl")
ensembl
## Object of class 'Mart':
## Using the ensembl BioMart database
## Using the dataset
Next you need to decide on a dataset. This can also be specified in the mart object that is created when you call the the useMart() method.
head(listDatasets(ensembl))
## dataset
## 1 oanatinus_gene_ensembl
## 2 cporcellus_gene_ensembl
## 3 gaculeatus_gene_ensembl
## 4 lafricana_gene_ensembl
## 5 itridecemlineatus_gene_ensembl
## 6 choffmanni_gene_ensembl
## description version
## 1 Ornithorhynchus anatinus genes (OANA5) OANA5
## 2 Cavia porcellus genes (cavPor3) cavPor3
## 3 Gasterosteus aculeatus genes (BROADS1) BROADS1
## 4 Loxodonta africana genes (loxAfr3) loxAfr3
## 5 Ictidomys tridecemlineatus genes (spetri2) spetri2
## 6 Choloepus hoffmanni genes (choHof1) choHof1
ensembl <- useMart("ensembl",dataset="hsapiens_gene_ensembl")
ensembl
## Object of class 'Mart':
## Using the ensembl BioMart database
## Using the hsapiens_gene_ensembl dataset
Next we need to think about filters and values. In the biomaRt package, filters are things that can be used with values to restrict or choose what comes back. So you might choose a filter of “affy_hg_u133_plus_2” to go with specific values. For example you might choose c(“202763_at”,“209310_s_at”,“207500_at”) to go with the filter “affy_hg_u133_plus_2”. Together these two things would request things that matched those probeset IDs on the platform listed as the filter. There is an accessor for the kinds of filters that are available from a given mart/dataset:
head(listFilters(ensembl))
## name description
## 1 chromosome_name Chromosome name
## 2 start Gene Start (bp)
## 3 end Gene End (bp)
## 4 band_start Band Start
## 5 band_end Band End
## 6 marker_start Marker Start
Also, you need to know about attributes. Attributes here mean the things that you want returned. So if you want to know the gene symbol or something like that. You would list that as an attribute. There are accessors to list the kinds of attributes you can look up too:
head(listAttributes(ensembl))
## name description
## 1 ensembl_gene_id Ensembl Gene ID
## 2 ensembl_transcript_id Ensembl Transcript ID
## 3 ensembl_peptide_id Ensembl Protein ID
## 4 ensembl_exon_id Ensembl Exon ID
## 5 description Description
## 6 chromosome_name Chromosome Name
Once you are done exploring and know what you want to extract, you can call the getBM method to get your data like this:
affyids=c("202763_at","209310_s_at","207500_at")
getBM(attributes=c('affy_hg_u133_plus_2', 'entrezgene'),
filters = 'affy_hg_u133_plus_2',
values = affyids, mart = ensembl)
## affy_hg_u133_plus_2 entrezgene
## 1 209310_s_at 837
## 2 207500_at 838
## 3 202763_at 836
Now what would you do if you didn't know what the possible values are for a given filter? Well you could just request all the possible values by not specifying the filter, and instead only specifying it as an attribute like this:
head(getBM(attributes='affy_hg_u133_plus_2', mart = ensembl))
## affy_hg_u133_plus_2
## 1 225862_at
## 2 1560940_at
## 3 1560941_a_at
## 4 240219_at
## 5 237885_at
## 6 222277_at
Of course if you find the standard biomaRt methods difficult to work with, you can now also use the standard select methods here.
[ Back to top ]
Exercise 11: Pull down GO terms for entrez gene id “1” from human by using the ensembl “hsapiens_gene_ensembl” dataset.
Exercise 12: Now compare the GO terms you just pulled down to the same GO terms from the org.Hs.eg.db package (which you can now retrieve using select()). What differences do you notice? Why do you suspect that is?
[ Back to top ]
There are many BSgenome packages in the repository too. These packages contain sequence data for sequenced organisms. You can load one of these packages just like this:
library(BSgenome.Hsapiens.UCSC.hg19)
ls(2)
## [1] "NP2009code" "attributePages" "columns" "exportFASTA"
## [5] "filterOptions" "filterType" "getBM" "getBMlist"
## [9] "getGene" "getLDS" "getSequence" "getXML"
## [13] "keys" "keytypes" "listAttributes" "listDatasets"
## [17] "listFilters" "listMarts" "select" "show"
## [21] "useDataset" "useMart"
Hsapiens
## Human genome
## |
## | organism: Homo sapiens (Human)
## | provider: UCSC
## | provider version: hg19
## | release date: Feb. 2009
## | release name: Genome Reference Consortium GRCh37
## |
## | single sequences (see '?seqnames'):
## | chr1 chr2 chr3
## | chr4 chr5 chr6
## | chr7 chr8 chr9
## | chr10 chr11 chr12
## | chr13 chr14 chr15
## | chr16 chr17 chr18
## | chr19 chr20 chr21
## | chr22 chrX chrY
## | chrM chr1_gl000191_random chr1_gl000192_random
## | chr4_ctg9_hap1 chr4_gl000193_random chr4_gl000194_random
## | chr6_apd_hap1 chr6_cox_hap2 chr6_dbb_hap3
## | chr6_mann_hap4 chr6_mcf_hap5 chr6_qbl_hap6
## | chr6_ssto_hap7 chr7_gl000195_random chr8_gl000196_random
## | chr8_gl000197_random chr9_gl000198_random chr9_gl000199_random
## | chr9_gl000200_random chr9_gl000201_random chr11_gl000202_random
## | chr17_ctg5_hap1 chr17_gl000203_random chr17_gl000204_random
## | chr17_gl000205_random chr17_gl000206_random chr18_gl000207_random
## | chr19_gl000208_random chr19_gl000209_random chr21_gl000210_random
## | chrUn_gl000211 chrUn_gl000212 chrUn_gl000213
## | chrUn_gl000214 chrUn_gl000215 chrUn_gl000216
## | chrUn_gl000217 chrUn_gl000218 chrUn_gl000219
## | chrUn_gl000220 chrUn_gl000221 chrUn_gl000222
## | chrUn_gl000223 chrUn_gl000224 chrUn_gl000225
## | chrUn_gl000226 chrUn_gl000227 chrUn_gl000228
## | chrUn_gl000229 chrUn_gl000230 chrUn_gl000231
## | chrUn_gl000232 chrUn_gl000233 chrUn_gl000234
## | chrUn_gl000235 chrUn_gl000236 chrUn_gl000237
## | chrUn_gl000238 chrUn_gl000239 chrUn_gl000240
## | chrUn_gl000241 chrUn_gl000242 chrUn_gl000243
## | chrUn_gl000244 chrUn_gl000245 chrUn_gl000246
## | chrUn_gl000247 chrUn_gl000248 chrUn_gl000249
## |
## | multiple sequences (see '?mseqnames'):
## | upstream1000 upstream2000 upstream5000
## |
## | (use the '$' or '[[' operator to access a given sequence)
The getSeq method is useful for extracting data from these pacakges. This method takes several arguments but the important ones are the 1st two. The 1st argument specifies the BSgenome object to use and the second argument (names) specifies what data you want back out. So for example, if you call it and give a character vector that names the seqnames for the object then you will get the sequences from those chromosomes as a DNAStringSet object.
seqNms <- seqnames(Hsapiens)
head(seqNms)
## [1] "chr1" "chr2" "chr3" "chr4" "chr5" "chr6"
getSeq(Hsapiens, seqNms[1:2])
## A DNAStringSet instance of length 2
## width seq
## [1] 249250621 NNNNNNNNNNNNNNNNNNNNNNNNNNNNN...NNNNNNNNNNNNNNNNNNNNNNNNNNNNN
## [2] 243199373 NNNNNNNNNNNNNNNNNNNNNNNNNNNNN...NNNNNNNNNNNNNNNNNNNNNNNNNNNNN
Whereas if you give the a GRanges object for the 2nd argument, you can instead get a DNA StringSet that corresponds to those ranges.
rngs <- GRanges(seqnames = c('chr1', 'chr4'), strand=c('+','-'),
ranges = IRanges(start=c(100000,300000),
end=c(100023,300037)))
rngs
## GRanges with 2 ranges and 0 metadata columns:
## seqnames ranges strand
## <Rle> <IRanges> <Rle>
## [1] chr1 [100000, 100023] +
## [2] chr4 [300000, 300037] -
## ---
## seqlengths:
## chr1 chr4
## NA NA
res <- getSeq(Hsapiens, rngs)
res
## A DNAStringSet instance of length 2
## width seq
## [1] 24 CACTAAGCACACAGAGAATAATGT
## [2] 38 GCTGGTCCCTTACTTCCAGTAGAAAAGACGTGTTCAGG
This can be a very powerful way to quickly get sequences of interest. And for more useful tools the BSgenome package also has useful functions for finding a pattern in a string set etc.
[ Back to top ]
Follow installation instructions to start using these packages. To install the annotations associated with the Affymetrix Human Genome U95 V 2.0, and with Gene Ontology, use
source("http://bioconductor.org/biocLite.R")
biocLite(c("hgu95av2.db", "GO.db"))
Package installation is required only once per R installation. View a full list of available software and annotation packages.
To use the AnnotationDbi
and GO.db
package, evaluate the commands
library(AnnotationDbi)
library(GO.db)
These commands are required once in each R session.
[ Back to top ]
Packages have extensive help pages, and include vignettes highlighting common use cases. The help pages and vignettes are available from within R. After loading a package, use syntax like
help(package="GO.db")
?select
to obtain an overview of help on the GO.db
package, and the select
method. The AnnotationDbi
package is used by most .db
packages. View the vignettes in the AnnotationDbi
package with
browseVignettes(package="AnnotationDbi")
To view vignettes (providing a more comprehensive introduction to
package functionality) in the AnnotationDbi
package. Use
help.start()
To open a web page containing comprehensive help resources.
[ Back to top ]
The following guides the user through key annotation packages. Users
interested in how to create custom chip packages should see the
vignettes in the AnnotationForge
package. There is additional
information in the AnnotationDbi
, OrganismDbi
and
GenomicFeatures
packages for how to use some of the extra tools
provided. You can also refer to the complete list of annotation
packages.
AnnotationDbi
package. This
package will be automatically installed for you if you install
another “.db” annotation package using biocLite(). It contains the
code to allow annotation mapping objects to be made and manipulated
as well as code to use the select methods etc..AnnotationDbi
package. These packages must be
upgraded before you attempt to update your custom chip packages as
they contain the source databases needed by the SQLForge code.[ Back to top ]
keys <- "MSX2"
columns <- c("ENTREZID","PROBEID", "CHR")
select(hgu95av2.db, keys, columns, keytype="SYMBOL")
## Warning: 'select' resulted in 1:many mapping between keys and return rows
## SYMBOL ENTREZID PROBEID CHR
## 1 MSX2 4488 40733_f_at 5
## 2 MSX2 4488 40734_r_at 5
Initially you might expect that hgu95av2.db will have less information in it. After all, it's an old Affymetrix platform that was developed before we even had a very complete human genome. So you might try something like this:
chipSymbols <- keys(hgu95av2.db, keytype="SYMBOL")
orgSymbols <- keys(org.Hs.eg.db, keytype="SYMBOL")
length(orgSymbols)
## [1] 47912
length(chipSymbols)
## [1] 47912
And you might feel confused and so you might try this:
dim(select(org.Hs.eg.db,orgSymbols, "ENTREZID", "SYMBOL"))
## Warning: 'select' resulted in 1:many mapping between keys and return rows
## [1] 47938 2
dim(select(hgu95av2.db,chipSymbols, "ENTREZID", "SYMBOL"))
## Warning: 'select' resulted in 1:many mapping between keys and return rows
## [1] 47938 2
And you might also have noticed this:
length(columns(org.Hs.eg.db)) < length(columns(hgu95av2.db))
## [1] TRUE
Well the answer you have in front of you is actually correct. There actually is more information available in the hgu95av2.db object than in the org.Hs.eg.db object. This is because even though the hgu95av2.db object technically can only have probes for some genes in the genome, it still (behind the scenes) retrieves data about gene names etc. from the org.Hs.eg.db package. So it effectively has access to all the data from the org package PLUS the probes for that platform and what those map to. So that means that for there will be information about many gene symbols that don't actually match up to any probeset Ids. And that is what we see if we use gene symbols to look up the probes Ids.
head(select(hgu95av2.db,chipSymbols, "PROBEID", "SYMBOL"))
## SYMBOL PROBEID
## 1 A1BG <NA>
## 2 A2M <NA>
## 3 A2MP1 <NA>
## 4 NAT1 38187_at
## 5 NAT2 38912_at
## 6 NATP <NA>
egr <- select(org.Hs.eg.db, orgSymbols, "ENTREZID", "SYMBOL")
## Warning: 'select' resulted in 1:many mapping between keys and return rows
length(egr$ENTREZID)
## [1] 47938
length(unique(egr$ENTREZID))
## [1] 47938
## VS:
length(egr$SYMBOL)
## [1] 47938
length(unique(egr$SYMBOL))
## [1] 47912
## So lets trap these symbols that are redundant and look more closely...
redund <- egr$SYMBOL
badSymbols <- redund[duplicated(redund)]
select(org.Hs.eg.db, badSymbols, "ENTREZID", "SYMBOL")
## Warning: 'select' resulted in 1:many mapping between keys and return rows
## SYMBOL ENTREZID
## 1 CSNK1E 1454
## 2 CSNK1E 102800317
## 3 HBD 3045
## 4 HBD 100187828
## 5 RNR1 4549
## 6 RNR1 6052
## 7 RNR2 4550
## 8 RNR2 6053
## 9 TEC 7006
## 10 TEC 100124696
## 11 MEMO1 7795
## 12 MEMO1 51072
## 13 DUX4 22947
## 14 DUX4 100653046
## 15 KIR3DL3 115653
## 16 KIR3DL3 100133046
## 17 MMD2 221938
## 18 MMD2 100505381
## 19 RP5-1043L13.1 284757
## 20 RP5-1043L13.1 729296
## 21 RP11-344E13.3 339260
## 22 RP11-344E13.3 440416
## 23 MIR642A 693227
## 24 MIR642A 102466336
## 25 LINC00623 728855
## 26 LINC00623 101929362
## 27 LINC00684 100129407
## 28 LINC00684 100132304
## 29 RP6-206I17.2 100130000
## 30 RP6-206I17.2 101929438
## 31 RP4-669L17.10 100132062
## 32 RP4-669L17.10 100132287
## 33 AC159540.1 100506076
## 34 AC159540.1 100506123
## 35 LSAMP-AS1 100506708
## 36 LSAMP-AS1 101926903
## 37 RP11-696N14.1 100507053
## 38 RP11-696N14.1 101929271
## 39 RP11-513O17.2 100507651
## 40 RP11-513O17.2 101929436
## 41 RP11-119F19.2 101060691
## 42 RP11-119F19.2 101929525
## 43 RP11-561O23.5 101926991
## 44 RP11-561O23.5 101929800
## 45 AC004893.11 101927550
## 46 AC004893.11 102723992
## 47 RP11-475O6.1 101927560
## 48 RP11-475O6.1 101927587
## 49 RP1-102G20.5 101928696
## 50 RP1-102G20.5 102724578
## 51 RP11-321G12.1 102723344
## 52 RP11-321G12.1 102723355
So to retrieve this information using select you need to do it like this:
res1 <- select(TxDb.Hsapiens.UCSC.hg19.knownGene,
keys(TxDb.Hsapiens.UCSC.hg19.knownGene, keytype="TXID"),
columns=c("GENEID","TXNAME","TXCHROM"), keytype="TXID")
head(res1)
## TXID GENEID TXNAME TXCHROM
## 1 1 100287102 uc001aaa.3 chr1
## 2 2 100287102 uc010nxq.1 chr1
## 3 3 100287102 uc010nxr.1 chr1
## 4 4 79501 uc001aal.1 chr1
## 5 5 <NA> uc001aaq.2 chr1
## 6 6 <NA> uc001aar.2 chr1
And to do it using transcripts you do it like this:
res2 <- transcripts(TxDb.Hsapiens.UCSC.hg19.knownGene,
columns = c("gene_id","tx_name"))
head(res2)
## GRanges with 6 ranges and 2 metadata columns:
## seqnames ranges strand | gene_id tx_name
## <Rle> <IRanges> <Rle> | <CharacterList> <character>
## [1] chr1 [ 11874, 14409] + | 100287102 uc001aaa.3
## [2] chr1 [ 11874, 14409] + | 100287102 uc010nxq.1
## [3] chr1 [ 11874, 14409] + | 100287102 uc010nxr.1
## [4] chr1 [ 69091, 70008] + | 79501 uc001aal.1
## [5] chr1 [321084, 321115] + | uc001aaq.2
## [6] chr1 [321146, 321207] + | uc001aar.2
## ---
## seqlengths:
## chr1 chr2 ... chrUn_gl000249
## 249250621 243199373 ... 38502
Notice that in the 2nd case we don't have to ask for the chromosome, as transcripts() returns a GRanges object, so the chromosome will automatically be returned as part of the object.
library(TxDb.Athaliana.BioMart.plantsmart16)
res <- transcripts(TxDb.Athaliana.BioMart.plantsmart16, columns = c("gene_id"))
You will notice that the gene ids for this package are TAIR locus IDs and are NOT entrez gene IDs like what you saw in the TxDb.Hsapiens.UCSC.hg19.knownGene package. It's important to always pay attention to the kind of gene id is being used by the TranscriptDb you are looking at.
library(Homo.sapiens)
keys <- keys(Homo.sapiens, keytype="TXID")
res1 <- select(Homo.sapiens,
keys= keys,
columns=c("SYMBOL","TXSTART","TXCHROM"), keytype="TXID")
head(res1)
And to do it using transcripts you do it like this:
library(Homo.sapiens)
res2 <- transcripts(Homo.sapiens, columns="SYMBOL")
head(res2)
## GRanges with 6 ranges and 1 metadata column:
## seqnames ranges strand | SYMBOL
## <Rle> <IRanges> <Rle> | <CharacterList>
## [1] chr1 [ 11874, 14409] + | DDX11L1
## [2] chr1 [ 11874, 14409] + | DDX11L1
## [3] chr1 [ 11874, 14409] + | DDX11L1
## [4] chr1 [ 69091, 70008] + | OR4F5
## [5] chr1 [321084, 321115] + | NA
## [6] chr1 [321146, 321207] + | NA
## ---
## seqlengths:
## chr1 chr2 ... chrUn_gl000249
## 249250621 243199373 ... 38502
columns(Homo.sapiens)
## [1] "GOID" "TERM" "ONTOLOGY" "DEFINITION"
## [5] "ENTREZID" "PFAM" "IPI" "PROSITE"
## [9] "ACCNUM" "ALIAS" "CHR" "CHRLOC"
## [13] "CHRLOCEND" "ENZYME" "MAP" "PATH"
## [17] "PMID" "REFSEQ" "SYMBOL" "UNIGENE"
## [21] "ENSEMBL" "ENSEMBLPROT" "ENSEMBLTRANS" "GENENAME"
## [25] "UNIPROT" "GO" "EVIDENCE" "GOALL"
## [29] "EVIDENCEALL" "ONTOLOGYALL" "OMIM" "UCSCKG"
## [33] "CDSID" "CDSNAME" "CDSCHROM" "CDSSTRAND"
## [37] "CDSSTART" "CDSEND" "EXONID" "EXONNAME"
## [41] "EXONCHROM" "EXONSTRAND" "EXONSTART" "EXONEND"
## [45] "GENEID" "TXID" "EXONRANK" "TXNAME"
## [49] "TXCHROM" "TXSTRAND" "TXSTART" "TXEND"
columns(org.Hs.eg.db)
## [1] "ENTREZID" "PFAM" "IPI" "PROSITE"
## [5] "ACCNUM" "ALIAS" "CHR" "CHRLOC"
## [9] "CHRLOCEND" "ENZYME" "MAP" "PATH"
## [13] "PMID" "REFSEQ" "SYMBOL" "UNIGENE"
## [17] "ENSEMBL" "ENSEMBLPROT" "ENSEMBLTRANS" "GENENAME"
## [21] "UNIPROT" "GO" "EVIDENCE" "ONTOLOGY"
## [25] "GOALL" "EVIDENCEALL" "ONTOLOGYALL" "OMIM"
## [29] "UCSCKG"
columns(TxDb.Hsapiens.UCSC.hg19.knownGene)
## [1] "CDSID" "CDSNAME" "CDSCHROM" "CDSSTRAND" "CDSSTART"
## [6] "CDSEND" "EXONID" "EXONNAME" "EXONCHROM" "EXONSTRAND"
## [11] "EXONSTART" "EXONEND" "GENEID" "TXID" "EXONRANK"
## [16] "TXNAME" "TXCHROM" "TXSTRAND" "TXSTART" "TXEND"
## You might also want to look at this:
transcripts(Homo.sapiens, columns=c("SYMBOL","CHRLOC"))
## Warning: 'select' resulted in 1:many mapping between keys and return rows
## Warning: 'select' resulted in 1:many mapping between keys and return rows
## GRanges with 82960 ranges and 3 metadata columns:
## seqnames ranges strand | CHRLOC
## <Rle> <IRanges> <Rle> | <IntegerList>
## [1] chr1 [ 11874, 14409] + | 11874
## [2] chr1 [ 11874, 14409] + | 11874
## [3] chr1 [ 11874, 14409] + | 11874
## [4] chr1 [ 69091, 70008] + | 69091
## [5] chr1 [321084, 321115] + | NA
## ... ... ... ... ... ...
## [82956] chrY [27605645, 27605678] - | NA
## [82957] chrY [27606394, 27606421] - | NA
## [82958] chrY [27607404, 27607432] - | NA
## [82959] chrY [27635919, 27635954] - | NA
## [82960] chrY [59358329, 59360854] - | NA
## CHRLOCCHR SYMBOL
## <CharacterList> <CharacterList>
## [1] 1 DDX11L1
## [2] 1 DDX11L1
## [3] 1 DDX11L1
## [4] 1 OR4F5
## [5] NA NA
## ... ... ...
## [82956] NA NA
## [82957] NA NA
## [82958] NA NA
## [82959] NA NA
## [82960] NA NA
## ---
## seqlengths:
## chr1 chr2 ... chrUn_gl000249
## 249250621 243199373 ... 38502
The key difference is that the TXSTART refers to the start of a transcript and originates in the TranscriptDb object from the TxDb.Hsapiens.UCSC.hg19.knownGene package, while the CHRLOC refers to the same thing but originates in the OrgDb object from the org.Hs.eg.db package. The point of origin is significant because the TranscriptDb object represents a transcriptome from UCSC and the OrgDb is primarily gene centric data that originates at NCBI. The upshot is that CHRLOC will not have as many regions represented as TXSTART, since there has to be an official gene for there to even be a record. The CHRLOC data is also locked in for org.Hs.eg.db as data for hg19, whereas you can swap in a different TranscriptDb object to match the genome you are using to make it hg18 etc. For these reasons, we strongly recommend using TXSTART instead of CHRLOC. Howeverm CHRLOC still remains in the org packages for historical reasons.
To find the keys that match, make use of the pattern and column arguments.
library(Homo.sapiens)
xk = head(keys(Homo.sapiens, keytype="ENTREZID", pattern="X", column="SYMBOL"))
xk
## [1] "100033409" "100033411" "100036519" "100038246" "100048904" "100048923"
select verifies the results
select(Homo.sapiens, xk, "SYMBOL", "ENTREZID")
## ENTREZID SYMBOL
## 1 100033409 OTX2P1
## 2 100033411 DUXB
## 3 100036519 FOXD4L2
## 4 100038246 TLX1NB
## 5 100048904 DDX39BP1
## 6 100048923 DDX39BP2
The 1st thing you need to do is look at the keytypes:
keytypes(ah)
## [1] "BiocVersion" "DataProvider" "Title"
## [4] "SourceFile" "Species" "SourceUrl"
## [7] "SourceVersion" "TaxonomyId" "Genome"
## [10] "Description" "Tags" "RDataClass"
## [13] "RDataPath" "Coordinate_1_based" "Maintainer"
## [16] "RDataVersion" "RDataDateAdded" "Recipe"
Then you want to look at possible values for DataProvider and for Genome.
keys(ah, keytype="DataProvider")
## [1] "EncodeDCC"
## [2] "ftp.ensembl.org"
## [3] "ftp://ftp.ncbi.nih.gov/snp"
## [4] "HAEMCODE"
## [5] "hgdownload.cse.ucsc.edu"
## [6] "http://inparanoid.sbc.su.se/download/current/Orthologs"
## [7] "RefNet"
head(keys(ah, keytype="Genome"))
## [1] "ailMel1" "anoCar1" "anoCar2" "AnoCar2.0" "anoGam1" "apiMel1"
filters(ah) <- NULL
filters(ah) <- list(Species="Homo sapiens",
DataProvider="hgdownload.cse.ucsc.edu",
Genome="hg19")
length(ah)
## [1] 118
This pulls down the oreganno annotations. Which are described on the UCSC site thusly: “This track displays literature-curated regulatory regions, transcription factor binding sites, and regulatory polymorphisms from ORegAnno (Open Regulatory Annotation). For more detailed information on a particular regulatory element, follow the link to ORegAnno from the details page.”
res <- ah$goldenpath.hg19.database.oreganno_0.0.1.RData
library("biomaRt")
ensembl <- useMart("ensembl",dataset="hsapiens_gene_ensembl")
ids=c("1")
getBM(attributes=c('go_id', 'entrezgene'),
filters = 'entrezgene',
values = ids, mart = ensembl)
## go_id entrezgene
## 1 GO:0008150 1
## 2 GO:0005576 1
## 3 GO:0005515 1
## 4 GO:0003674 1
library(org.Hs.eg.db)
ids=c("1")
select(org.Hs.eg.db, keys=ids, columns="GO", keytype="ENTREZID")
## Warning: 'select' resulted in 1:many mapping between keys and return rows
## ENTREZID GO EVIDENCE ONTOLOGY
## 1 1 GO:0003674 ND MF
## 2 1 GO:0005576 IDA CC
## 3 1 GO:0008150 ND BP
## 4 1 GO:0070062 IDA CC
## 5 1 GO:0072562 IDA CC
When this exercise was written, there was a different number of GO terms returned from biomaRt than from org.Hs.eg.db. This may not always be true in the future though as both of these resources are updated. It is expected however that this web service, (which is updated continuously) will fall in and out of sync with the org.Hs.eg.db package (which is updated twice a year). This is an important difference as each approach has different advantages and disadvantages. The advantage to updating continuously is that you always have the very latest annotations which are frequently different for something like GO terms. The advantage to using a package is that the results are frozen to a release of Bioconductor. And this can help you to get the same answers that you get today (reproducibility), a few years from now.
[ Back to top ]
sessionInfo()
## R version 3.1.0 (2014-04-10)
## Platform: x86_64-apple-darwin10.8.0 (64-bit)
##
## locale:
## [1] C
##
## attached base packages:
## [1] parallel stats graphics grDevices utils datasets methods
## [8] base
##
## other attached packages:
## [1] TxDb.Athaliana.BioMart.plantsmart16_2.9.0
## [2] biomaRt_2.20.0
## [3] Homo.sapiens_1.1.2
## [4] OrganismDbi_1.6.0
## [5] GO.db_2.14.0
## [6] hgu95av2.db_2.14.0
## [7] ensemblVEP_1.4.0
## [8] BSgenome.Hsapiens.UCSC.hg19_1.3.99
## [9] BSgenome_1.32.0
## [10] org.Mm.eg.db_2.14.0
## [11] org.Hs.eg.db_2.14.0
## [12] RSQLite_0.11.4
## [13] DBI_0.2-7
## [14] TxDb.Mmusculus.UCSC.mm10.ensGene_2.14.0
## [15] TxDb.Hsapiens.UCSC.hg19.knownGene_2.14.0
## [16] GenomicFeatures_1.16.2
## [17] AnnotationDbi_1.26.0
## [18] Biobase_2.24.0
## [19] AnnotationHub_1.4.0
## [20] VariantAnnotation_1.10.5
## [21] Rsamtools_1.16.1
## [22] Biostrings_2.32.1
## [23] XVector_0.4.0
## [24] GenomicRanges_1.16.3
## [25] GenomeInfoDb_1.0.2
## [26] IRanges_1.22.9
## [27] BiocGenerics_0.10.0
##
## loaded via a namespace (and not attached):
## [1] BBmisc_1.7 BatchJobs_1.2
## [3] BiocInstaller_1.14.2 BiocParallel_0.6.1
## [5] Category_2.30.0 GSEABase_1.26.0
## [7] GenomicAlignments_1.0.2 MASS_7.3-33
## [9] Matrix_1.1-3 RBGL_1.40.0
## [11] RColorBrewer_1.0-5 RCurl_1.95-4.1
## [13] RJSONIO_1.2-0.2 Rcpp_0.11.2
## [15] XML_3.98-1.1 annotate_1.42.0
## [17] bitops_1.0-6 brew_1.0-6
## [19] caTools_1.17 checkmate_1.1
## [21] codetools_0.2-8 colorspace_1.2-4
## [23] digest_0.6.4 evaluate_0.5.5
## [25] fail_1.2 foreach_1.4.2
## [27] formatR_0.10 genefilter_1.46.1
## [29] ggplot2_1.0.0 graph_1.42.0
## [31] grid_3.1.0 gridSVG_1.4-0
## [33] gtable_0.1.2 htmltools_0.2.4
## [35] httpuv_1.3.0 httr_0.3
## [37] interactiveDisplay_1.2.0 iterators_1.0.7
## [39] knitr_1.6 lattice_0.20-29
## [41] markdown_0.7 munsell_0.4.2
## [43] plyr_1.8.1 proto_0.3-10
## [45] reshape2_1.4 rjson_0.2.14
## [47] rtracklayer_1.24.2 scales_0.2.4
## [49] sendmailR_1.1-2 shiny_0.10.0
## [51] splines_3.1.0 stats4_3.1.0
## [53] stringr_0.6.2 survival_2.37-7
## [55] tools_3.1.0 xtable_1.7-3
## [57] zlibbioc_1.10.0
[ Back to top ]