Package version: biomaRt 2.30.0

Contents

1 Introduction

In recent years a wealth of biological data has become available in public data repositories. Easy access to these valuable data resources and firm integration with data analysis is needed for comprehensive bioinformatics data analysis. The biomaRt package, provides an interface to a growing collection of databases implementing the BioMart software suite. The package enables retrieval of large amounts of data in a uniform way without the need to know the underlying database schemas or write complex SQL queries. Examples of BioMart databases are Ensembl, Uniprot and HapMap. These major databases give biomaRt users direct access to a diverse set of data and enable a wide range of powerful online queries from R.

2 Selecting a BioMart database and dataset

Every analysis with biomaRt starts with selecting a BioMart database to use. A first step is to check which BioMart web services are available. The function listMarts() will display all available BioMart web services

library("biomaRt")
listMarts()
##                biomart               version
## 1 ENSEMBL_MART_ENSEMBL      Ensembl Genes 86
## 2   ENSEMBL_MART_MOUSE      Mouse strains 86
## 3     ENSEMBL_MART_SNP  Ensembl Variation 86
## 4 ENSEMBL_MART_FUNCGEN Ensembl Regulation 86
## 5    ENSEMBL_MART_VEGA               Vega 66

Note: if the function useMart() runs into proxy problems you should set your proxy first before calling any biomaRt functions.
You can do this using the Sys.putenv command:

Sys.setenv("http_proxy" = "http://my.proxy.org:9999")

Some users have reported that the workaround above does not work, in this case an alternative proxy solution below can be tried:

options(RCurlOptions = list(proxy="uscache.kcc.com:80",proxyuserpwd="------:-------"))

The useMart() function can now be used to connect to a specified BioMart database, this must be a valid name given by listMarts(). In the next example we choose to query the Ensembl BioMart database.

ensembl=useMart("ensembl")

BioMart databases can contain several datasets, for Ensembl every species is a different dataset. In a next step we look at which datasets are available in the selected BioMart by using the function listDatasets().

listDatasets(ensembl)
##                           dataset                                description           version
## 1          oanatinus_gene_ensembl     Ornithorhynchus anatinus genes (OANA5)             OANA5
## 2         cporcellus_gene_ensembl            Cavia porcellus genes (cavPor3)           cavPor3
## 3         gaculeatus_gene_ensembl     Gasterosteus aculeatus genes (BROADS1)           BROADS1
## 4  itridecemlineatus_gene_ensembl Ictidomys tridecemlineatus genes (spetri2)           spetri2
## 5          lafricana_gene_ensembl         Loxodonta africana genes (loxAfr3)           loxAfr3
## 6         choffmanni_gene_ensembl        Choloepus hoffmanni genes (choHof1)           choHof1
## 7          csavignyi_gene_ensembl             Ciona savignyi genes (CSAV2.0)           CSAV2.0
## 8             fcatus_gene_ensembl        Felis catus genes (Felis_catus_6.2)   Felis_catus_6.2
## 9        rnorvegicus_gene_ensembl         Rattus norvegicus genes (Rnor_6.0)          Rnor_6.0
## 10         psinensis_gene_ensembl     Pelodiscus sinensis genes (PelSin_1.0)        PelSin_1.0
## 11          cjacchus_gene_ensembl  Callithrix jacchus genes (C_jacchus3.2.1)    C_jacchus3.2.1
## 12        ttruncatus_gene_ensembl         Tursiops truncatus genes (turTru1)           turTru1
## 13       scerevisiae_gene_ensembl   Saccharomyces cerevisiae genes (R64-1-1)           R64-1-1
## 14          celegans_gene_ensembl    Caenorhabditis elegans genes (WBcel235)          WBcel235
## 15          csabaeus_gene_ensembl      Chlorocebus sabaeus genes (ChlSab1.1)         ChlSab1.1
## 16        oniloticus_gene_ensembl    Oreochromis niloticus genes (Orenil1.0)         Orenil1.0
## 17        amexicanus_gene_ensembl       Astyanax mexicanus genes (AstMex102)         AstMex102
## 18         trubripes_gene_ensembl            Takifugu rubripes genes (FUGU4)             FUGU4
## 19          pmarinus_gene_ensembl    Petromyzon marinus genes (Pmarinus_7.0)      Pmarinus_7.0
## 20        eeuropaeus_gene_ensembl       Erinaceus europaeus genes (HEDGEHOG)          HEDGEHOG
## 21       falbicollis_gene_ensembl     Ficedula albicollis genes (FicAlb_1.4)        FicAlb_1.4
## 22         etelfairi_gene_ensembl           Echinops telfairi genes (TENREC)            TENREC
## 23     cintestinalis_gene_ensembl              Ciona intestinalis genes (KH)                KH
## 24      ptroglodytes_gene_ensembl         Pan troglodytes genes (CHIMP2.1.4)        CHIMP2.1.4
## 25       nleucogenys_gene_ensembl        Nomascus leucogenys genes (Nleu1.0)           Nleu1.0
## 26           sscrofa_gene_ensembl             Sus scrofa genes (Sscrofa10.2)       Sscrofa10.2
## 27        ocuniculus_gene_ensembl    Oryctolagus cuniculus genes (OryCun2.0)         OryCun2.0
## 28     dnovemcinctus_gene_ensembl     Dasypus novemcinctus genes (Dasnov3.0)         Dasnov3.0
## 29         pcapensis_gene_ensembl          Procavia capensis genes (proCap1)           proCap1
## 30          tguttata_gene_ensembl    Taeniopygia guttata genes (taeGut3.2.4)       taeGut3.2.4
## 31        mlucifugus_gene_ensembl         Myotis lucifugus genes (Myoluc2.0)         Myoluc2.0
## 32          hsapiens_gene_ensembl             Homo sapiens genes (GRCh38.p7)         GRCh38.p7
## 33          pformosa_gene_ensembl      Poecilia formosa genes (PoeFor_5.1.2)      PoeFor_5.1.2
## 34        tbelangeri_gene_ensembl         Tupaia belangeri genes (TREESHREW)         TREESHREW
## 35             mfuro_gene_ensembl Mustela putorius furo genes (MusPutFur1.0)      MusPutFur1.0
## 36           ggallus_gene_ensembl    Gallus gallus genes (Gallus_gallus-5.0) Gallus_gallus-5.0
## 37       xtropicalis_gene_ensembl         Xenopus tropicalis genes (JGI_4.2)           JGI_4.2
## 38         ecaballus_gene_ensembl             Equus caballus genes (EquCab2)           EquCab2
## 39           pabelii_gene_ensembl                 Pongo abelii genes (PPYG2)             PPYG2
## 40            drerio_gene_ensembl                 Danio rerio genes (GRCz10)            GRCz10
## 41        xmaculatus_gene_ensembl  Xiphophorus maculatus genes (Xipmac4.4.2)       Xipmac4.4.2
## 42     tnigroviridis_gene_ensembl  Tetraodon nigroviridis genes (TETRAODON8)        TETRAODON8
## 43        lchalumnae_gene_ensembl        Latimeria chalumnae genes (LatCha1)           LatCha1
## 44      amelanoleuca_gene_ensembl     Ailuropoda melanoleuca genes (ailMel1)           ailMel1
## 45          mmulatta_gene_ensembl          Macaca mulatta genes (Mmul_8.0.1)        Mmul_8.0.1
## 46         pvampyrus_gene_ensembl          Pteropus vampyrus genes (pteVam1)           pteVam1
## 47           panubis_gene_ensembl             Papio anubis genes (PapAnu2.0)         PapAnu2.0
## 48        mdomestica_gene_ensembl      Monodelphis domestica genes (BROADO5)           BROADO5
## 49     acarolinensis_gene_ensembl      Anolis carolinensis genes (AnoCar2.0)         AnoCar2.0
## 50            vpacos_gene_ensembl              Vicugna pacos genes (vicPac1)           vicPac1
## 51         tsyrichta_gene_ensembl           Tarsius syrichta genes (tarSyr1)           tarSyr1
## 52        ogarnettii_gene_ensembl         Otolemur garnettii genes (OtoGar3)           OtoGar3
## 53     dmelanogaster_gene_ensembl      Drosophila melanogaster genes (BDGP6)             BDGP6
## 54          mmurinus_gene_ensembl        Microcebus murinus genes (Mmur_2.0)          Mmur_2.0
## 55         loculatus_gene_ensembl       Lepisosteus oculatus genes (LepOcu1)           LepOcu1
## 56          olatipes_gene_ensembl            Oryzias latipes genes (MEDAKA1)           MEDAKA1
## 57         oprinceps_gene_ensembl             Ochotona princeps genes (pika)              pika
## 58          ggorilla_gene_ensembl          Gorilla gorilla genes (gorGor3.1)         gorGor3.1
## 59            dordii_gene_ensembl            Dipodomys ordii genes (dipOrd1)           dipOrd1
## 60            oaries_gene_ensembl                Ovis aries genes (Oar_v3.1)          Oar_v3.1
## 61         mmusculus_gene_ensembl             Mus musculus genes (GRCm38.p4)         GRCm38.p4
## 62        mgallopavo_gene_ensembl           Meleagris gallopavo genes (UMD2)              UMD2
## 63           gmorhua_gene_ensembl               Gadus morhua genes (gadMor1)           gadMor1
## 64          saraneus_gene_ensembl        Sorex araneus genes (COMMON_SHREW1)     COMMON_SHREW1
## 65    aplatyrhynchos_gene_ensembl    Anas platyrhynchos genes (BGI_duck_1.0)      BGI_duck_1.0
## 66         sharrisii_gene_ensembl      Sarcophilus harrisii genes (DEVIL7.0)          DEVIL7.0
## 67          meugenii_gene_ensembl          Macropus eugenii genes (Meug_1.0)          Meug_1.0
## 68           btaurus_gene_ensembl                  Bos taurus genes (UMD3.1)            UMD3.1
## 69       cfamiliaris_gene_ensembl         Canis familiaris genes (CanFam3.1)         CanFam3.1

To select a dataset we can update the Mart object using the function useDataset(). In the example below we choose to use the hsapiens dataset.

ensembl = useDataset("hsapiens_gene_ensembl",mart=ensembl)

Or alternatively if the dataset one wants to use is known in advance, we can select a BioMart database and dataset in one step by:

ensembl = useMart("ensembl",dataset="hsapiens_gene_ensembl")

3 How to build a biomaRt query

The getBM() function has three arguments that need to be introduced: filters, attributes and values. Filters define a restriction on the query. For example you want to restrict the output to all genes located on the human X chromosome then the filter chromosome_name can be used with value ‘X’. The listFilters() function shows you all available filters in the selected dataset.

filters = listFilters(ensembl)
filters[1:5,]
##              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

Attributes define the values we are interested in to retrieve. For example we want to retrieve the gene symbols or chromosomal coordinates. The listAttributes() function displays all available attributes in the selected dataset.

attributes = listAttributes(ensembl)
attributes[1:5,]
##                    name           description         page
## 1       ensembl_gene_id       Ensembl Gene ID feature_page
## 2 ensembl_transcript_id Ensembl Transcript ID feature_page
## 3    ensembl_peptide_id    Ensembl Protein ID feature_page
## 4       ensembl_exon_id       Ensembl Exon ID feature_page
## 5           description           Description feature_page

The getBM() function is the main query function in biomaRt. It has four main arguments:

Note: for some frequently used queries to Ensembl, wrapper functions are available: getGene() and getSequence(). These functions call the getBM() function with hard coded filter and attribute names.

Now that we selected a BioMart database and dataset, and know about attributes, filters, and the values for filters; we can build a biomaRt query. Let’s make an easy query for the following problem: We have a list of Affymetrix identifiers from the u133plus2 platform and we want to retrieve the corresponding EntrezGene identifiers using the Ensembl mappings.

The u133plus2 platform will be the filter for this query and as values for this filter we use our list of Affymetrix identifiers. As output (attributes) for the query we want to retrieve the EntrezGene and u133plus2 identifiers so we get a mapping of these two identifiers as a result. The exact names that we will have to use to specify the attributes and filters can be retrieved with the listAttributes() and listFilters() function respectively. Let’s now run the query:

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           207500_at        838
## 2           202763_at        836
## 3         209310_s_at        837

4 Examples of biomaRt queries

In the sections below a variety of example queries are described. Every example is written as a task, and we have to come up with a biomaRt solution to the problem.

4.1 Annotate a set of Affymetrix identifiers with HUGO symbol and chromosomal locations of corresponding genes

We have a list of Affymetrix hgu133plus2 identifiers and we would like to retrieve the HUGO gene symbols, chromosome names, start and end positions and the bands of the corresponding genes. The listAttributes() and the listFilters() functions give us an overview of the available attributes and filters and we look in those lists to find the corresponding attribute and filter names we need. For this query we’ll need the following attributes: hgnc_symbol, chromsome_name, start_position, end_position, band and affy_hg_u133_plus_2 (as we want these in the output to provide a mapping with our original Affymetrix input identifiers. There is one filter in this query which is the affy_hg_u133_plus_2 filter as we use a list of Affymetrix identifiers as input. Putting this all together in the getBM() and performing the query gives:

affyids=c("202763_at","209310_s_at","207500_at")
getBM(attributes = c('affy_hg_u133_plus_2', 'hgnc_symbol', 'chromosome_name',
                   'start_position', 'end_position', 'band'),
      filters = 'affy_hg_u133_plus_2', 
      values = affyids, 
      mart = ensembl)
##   affy_hg_u133_plus_2 hgnc_symbol chromosome_name start_position end_position  band
## 1           207500_at       CASP5              11      104994235    105023168 q22.3
## 2           202763_at       CASP3               4      184627696    184649509 q35.1
## 3         209310_s_at       CASP4              11      104942866    104969436 q22.3

4.2 Annotate a set of EntrezGene identifiers with GO annotation

In this task we start out with a list of EntrezGene identiers and we want to retrieve GO identifiers related to biological processes that are associated with these entrezgene identifiers. Again we look at the output of listAttributes() and listFilters() to find the filter and attributes we need. Then we construct the following query:

entrez=c("673","837")
goids = getBM(attributes = c('entrezgene', 'go_id'), 
              filters = 'entrezgene', 
              values = entrez, 
              mart = ensembl)
head(goids)
##   entrezgene      go_id
## 1        673           
## 2        673 GO:0015758
## 3        673 GO:0071277
## 4        673 GO:0090150
## 5        673 GO:0009887
## 6        673 GO:0070374

4.3 Retrieve all HUGO gene symbols of genes that are located on chromosomes 17,20 or Y, and are associated with specific GO terms

The GO terms we are interested in are: GO:0051330, GO:0000080, GO:0000114, GO:0000082. The key to performing this query is to understand that the getBM() function enables you to use more than one filter at the same time. In order to do this, the filter argument should be a vector with the filter names. The values should be a list, where the first element of the list corresponds to the first filter and the second list element to the second filter and so on. The elements of this list are vectors containing the possible values for the corresponding filters.

 go=c("GO:0051330","GO:0000080","GO:0000114","GO:0000082")
 chrom=c(17,20,"Y")
 getBM(attributes= "hgnc_symbol",
        filters=c("go_id","chromosome_name"),
        values=list(go, chrom), mart=ensembl)
##   hgnc_symbol
## 1        CDC6
## 2        CDK3
## 3        PCNA
## 4       CRLF3
## 5        RPA1
## 6        MCM8
## 7     RPS6KB1

4.4 Annotate set of idenfiers with INTERPRO protein domain identifiers

In this example we want to annotate the following two RefSeq identifiers: NM_005359 and NM_000546 with INTERPRO protein domain identifiers and a description of the protein domains.

refseqids = c("NM_005359","NM_000546")
ipro = getBM(attributes=c("refseq_mrna","interpro","interpro_description"), 
             filters="refseq_mrna",
             values=refseqids, 
             mart=ensembl)
ipro
##    refseq_mrna  interpro                                   interpro_description
## 1    NM_000546 IPR002117                           p53 tumour suppressor family
## 2    NM_000546 IPR008967             p53-like transcription factor, DNA-binding
## 3    NM_000546 IPR010991                            p53, tetramerisation domain
## 4    NM_000546 IPR011615                                p53, DNA-binding domain
## 5    NM_000546 IPR012346 p53/RUNT-type transcription factor, DNA-binding domain
## 6    NM_000546 IPR013872                             p53 transactivation domain
## 7    NM_005359 IPR001132                              SMAD domain, Dwarfin-type
## 8    NM_005359 IPR003619                           MAD homology 1, Dwarfin-type
## 9    NM_005359 IPR008984                                        SMAD/FHA domain
## 10   NM_005359 IPR013019                                      MAD homology, MH1
## 11   NM_005359 IPR013790                                                Dwarfin
## 12   NM_005359 IPR017855                                       SMAD domain-like

4.5 Select all Affymetrix identifiers on the hgu133plus2 chip and Ensembl gene identifiers for genes located on chromosome 16 between basepair 1100000 and 1250000.

In this example we will again use multiple filters: chromosome_name, start, and end as we filter on these three conditions. Note that when a chromosome name, a start position and an end position are jointly used as filters, the BioMart webservice interprets this as return everything from the given chromosome between the given start and end positions.

getBM(attributes = c('affy_hg_u133_plus_2','ensembl_gene_id'), 
      filters = c('chromosome_name','start','end'),
      values = list(16,1100000,1250000), 
      mart = ensembl)
##    affy_hg_u133_plus_2 ensembl_gene_id
## 1                      ENSG00000260702
## 2            215502_at ENSG00000260532
## 3                      ENSG00000273551
## 4            205845_at ENSG00000196557
## 5                      ENSG00000196557
## 6                      ENSG00000260403
## 7                      ENSG00000259910
## 8                      ENSG00000261294
## 9          220339_s_at ENSG00000116176
## 10                     ENSG00000277010
## 11         217023_x_at ENSG00000197253
## 12         210084_x_at ENSG00000197253
## 13         215382_x_at ENSG00000197253
## 14         216474_x_at ENSG00000197253
## 15         207134_x_at ENSG00000197253
## 16         205683_x_at ENSG00000197253
## 17         217023_x_at ENSG00000172236
## 18         210084_x_at ENSG00000172236
## 19         215382_x_at ENSG00000172236
## 20         207741_x_at ENSG00000172236
## 21         216474_x_at ENSG00000172236
## 22         207134_x_at ENSG00000172236
## 23         205683_x_at ENSG00000172236

4.6 Retrieve all entrezgene identifiers and HUGO gene symbols of genes which have a “MAP kinase activity” GO term associated with it.

The GO identifier for MAP kinase activity is GO:0004707. In our query we will use go_id as our filter, and entrezgene and hgnc_symbol as attributes. Here’s the query:

getBM(attributes = c('entrezgene','hgnc_symbol'), 
      filters = 'go_id', 
      values = 'GO:0004707', 
      mart = ensembl)
##    entrezgene hgnc_symbol
## 1      225689      MAPK15
## 2       51701         NLK
## 3        5594       MAPK1
## 4        5597       MAPK6
## 5        6300      MAPK12
## 6        5600      MAPK11
## 7        5603      MAPK13
## 8        1432      MAPK14
## 9        5596       MAPK4
## 10       5598       MAPK7
## 11       5595       MAPK3

4.7 Given a set of EntrezGene identifiers, retrieve 100bp upstream promoter sequences

All sequence related queries to Ensembl are available through the getSequence() wrapper function. getBM() can also be used directly to retrieve sequences but this can get complicated so using getSequence is recommended.

Sequences can be retrieved using the getSequence() function either starting from chromosomal coordinates or identifiers.
The chromosome name can be specified using the chromosome argument. The start and end arguments are used to specify start and end positions on the chromosome. The type of sequence returned can be specified by the seqType argument which takes the following values:

In MySQL mode the getSequence() function is more limited and the sequence that is returned is the 5’ to 3’+ strand of the genomic sequence, given a chromosome, as start and an end position.

This task requires us to retrieve 100bp upstream promoter sequences from a set of EntrzGene identifiers. The type argument in getSequence() can be thought of as the filter in this query and uses the same input names given by listFilters(). In our query we use entrezgene for the type argument. Next we have to specify which type of sequences we want to retrieve, here we are interested in the sequences of the promoter region, starting right next to the coding start of the gene. Setting the seqType to coding_gene_flank will give us what we need. The upstream argument is used to specify how many bp of upstream sequence we want to retrieve, here we’ll retrieve a rather short sequence of 100bp. Putting this all together in getSequence() gives:

entrez=c("673","7157","837")
getSequence(id = entrez, 
            type="entrezgene",
            seqType="coding_gene_flank",
            upstream=100, 
            mart=ensembl) 
##                                                                                      coding_gene_flank entrezgene
## 1 CCTCCGCCTCCGCCTCCGCCTCCGCCTCCCCCAGCTCTCCGCCTCCCTTCCCCCTCCCCGCCCGACAGCGGCCGCTCGGGCCCCGGCTCTCGGTTATAAG        673
## 2 CACGTTTCCGCCCTTTGCAATAAGGAAATACATAGTTTACTTTCATTTTTGACTCTGAGGCTCTTTCCAACGCTGTAAAAAAGGACAGAGGCTGTTCCCT        837
## 3 TCCTTCTCTGCAGGCCCAGGTGACCCAGGGTTGGAAGTGTCTCATGCTGGATCCCCACTTTTCCTCTTGCAGCAGCCAGACTGCCTTCCGGGTCACTGCC       7157

4.8 Retrieve all 5’ UTR sequences of all genes that are located on chromosome 3 between the positions 185,514,033 and 185,535,839

As described in the provious task getSequence can also use chromosomal coordinates to retrieve sequences of all genes that lie in the given region. We also have to specify which type of identifier we want to retrieve together with the sequences, here we choose for entrezgene identifiers.

utr5 = getSequence(chromosome=3, start=185514033, end=185535839,
                   type="entrezgene",
                   seqType="5utr", 
                   mart=ensembl)
utr5
##                                                                                                                                             5utr
## 1                                                                                                                           Sequence unavailable
## 2                                                                                                        ATTCTTGTGAATGTGACACACGATCTCTCCAGTTTCCAT
## 3                                                        TGAGCAAAATCCCACAGTGGAAACTCTTAAGCCTCTGCGAAGTAAATCATTCTTGTGAATGTGACACACGATCTCTCCAGTTTCCAT
## 4 AGTCCCTAGGGAACTTCCTGTTGTCACCACACCTCTGAGTCGTCTGAGCTCACTGTGAGCAAAATCCCACAGTGGAAACTCTTAAGCCTCTGCGAAGTAAATCATTCTTGTGAATGTGACACACGATCTCTCCAGTTTCCAT
##   entrezgene
## 1     200879
## 2     200879
## 3     200879
## 4     200879

4.9 Retrieve protein sequences for a given list of EntrezGene identifiers

In this task the type argument specifies which type of identifiers we are using. To get an overview of other valid identifier types we refer to the listFilters() function.

protein = getSequence(id=c(100, 5728),
                      type="entrezgene",
                      seqType="peptide", 
                      mart=ensembl)
protein
##                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                              peptide
## 1                                                                                                                                                                                                                                                                                                                           MAQTPAFDKPKVELHVHLDGSIKPETILYYGRRRGIALPANTAEGLLNVIGMDKPLTLPDFLAKFDYYMPAIAGCREAIKRIAYEFVEMKAKEGVVYVEVRYSPHLLANSKVEPIPWNQAEGDLTPDEVVALVGQGLQEGERDFGVKARSILCCMRHQPNWSPKVVELCKKYQQQTVVAIDLAGDETIPGSSLLPGHVQAYQEAVKSGIHRTVHAGEVGSAEVVKEAVDILKTERLGHGYHTLEDQALYNRLRQENMHFEAQK*
## 2                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                       MAQTPAFDKPKVELHVHLDGSIKPETILYYGRRRGIALPANTAEGLLNVIGMDKPLTLPDFLAKFDYYMPAIARL*
## 3                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                               Sequence unavailable
## 4                                                                                                                                                                                                                                                                                                                                                                                                                                         ALLFHKMMFETIPMFSGGTCNPQFVVCQLKVKIYSSNSGPTRREDKFMYFEFPQPLPVCGDIKVEFFHKQNKMLKKDKMFHFWVNTFFIPGPEETSEKVENGSLCDQEIDSICSIERADNDKEYLVLTLTKNDLDKANKDKANRYFSPNFKVS*
## 5                                                                                                                                                                                                                       MAQTPAFDKPKVELHVHLDGSIKPETILYYGRRRGIALPANTAEGLLNVIGMDKPLTLPDFLAKFDYYMPAIAGCREAIKRIAYEFVEMKAKEGVVYVEVRYSPHLLANSKVEPIPWNQAEGDLTPDEVVALVGQGLQEGERDFGVKARSILCCMRHQPNWSPKVVELCKKYQQQTVVAIDLAGDETIPGSSLLPGHVQAYQEAVKSGIHRTVHAGEVGSAEVVKEAVDILKTERLGHGYHTLEDQALYNRLRQENMHFEICPWSSYLTGAWKPDTEHAVIRLKNDQANYSLNTDDPLIFKSTLDTDYQMTKRDMGFTEEEFKRLNINAAKSSFLPEDEKRELLDLLYKAYGMPPSASAGQNL*
## 6                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                                               Sequence unavailable
## 7  LERGGEAAAAAAAAAAAPGRGSESPVTISRAGNAGELVSPLLLPPTRRRRRRHIQGPGPVLNLPSAAAAPPVARAPEAAGGGSRSEDYSSSPHSAAAAARPLAAEEKQAQSLQPSSSRRSSHYPAAVQSQAAAERGASATAKSRAISILQKKPRHQQLLPSLSSFFFSHRLPDMTAIIKEIVSRNKRRYQEDGFDLDLTYIYPNIIAMGFPAERLEGVYRNNIDDVVRFLDSKHKNHYKIYNLCAERHYDTAKFNCRVAQYPFEDHNPPQLELIKPFCEDLDQWLSEDDNHVAAIHCKAGKGRTGVMICAYLLHRGKFLKAQEALDFYGEVRTRDKKGVTIPSQRRYVYYYSYLLKNHLDYRPVALLFHKMMFETIPMFSGGTCNPQFVVCQLKVKIYSSNSGPTRREDKFMYFEFPQPLPVCGDIKVEFFHKQNKMLKKDKMFHFWVNTFFIPGPEETSEKVENGSLCDQEIDSICSIERADNDKEYLVLTLTKNDLDKANKDKANRYFSPNFKVKLYFTKTVEEPSNPEASSSTSVTPDVSDNEPDHYRYSDTTDSDPENEPFDEDQHTQITKV*
## 8                                                                                                                                                                                                                                               MAQTPAFDKPKVELHVHLDGSIKPETILYYGRRRGIALPANTAEGLLNVIGMDKPLTLPDFLAKFDYYMPAIAGCREAIKRIAYEFVEMKAKEGVVYVEVRYSPHLLANSKVEPIPWNQAEGDLTPDEVVALVGQGLQEGERDFGVKARSILCCMRHQPNWSPKVVELCKKYQQQTVVAIDLAGDETIPGSSLLPGHVQAYQAVDILKTERLGHGYHTLEDQALYNRLRQENMHFEICPWSSYLTGAWKPDTEHAVIRLKNDQANYSLNTDDPLIFKSTLDTDYQMTKRDMGFTEEEFKRLNINAAKSSFLPEDEKRELLDLLYKAYGMPPSASAGQNL*
## 9                                                                                                                                                                               MTAIIKEIVSRNKRRYQEDGFDLDLTYIYPNIIAMGFPAERLEGVYRNNIDDVVRFLDSKHKNHYKIYNLCAERHYDTAKFNCRVAQYPFEDHNPPQLELIKPFCEDLDQWLSEDDNHVAAIHCKAGKGRTGVMICAYLLHRGKFLKAQEALDFYGEVRTRDKKGVTIPSQRRYVYYYSYLLKNHLDYRPVALLFHKMMFETIPMFSGGTCNPQFVVCQLKVKIYSSNSGPTRREDKFMYFEFPQPLPVCGDIKVEFFHKQNKMLKKDKMFHFWVNTFFIPGPEETSEKVENGSLCDQEIDSICSIERADNDKEYLVLTLTKNDLDKANKDKANRYFSPNFKVKLYFTKTVEEPSNPEASSSTSVTPDVSDNEPDHYRYSDTTDSDPENEPFDEDQHTQITKV*
## 10                                                                                                                                                                                                                                                                                                                                                                                                                                     MGFPAERLEGVYRNNIDDVVRFLDSKHKNHYKIYNLCAERHYDTAKFNCRVAQYPFEDHNPPQLELIKPFCEDLDQWLSEDDNHVAAIHCKAGKGRTGVMICAYLLHRGKFLKAQEALDFYGEVRTRDKKADPTGGIPDKGIIVIGDGSSMDVIAP*
## 11         RLGDSALAPRATAAAPGRGSESPVTISRAGNAGELVSPLLLPPTIQGPGPVLNLPSAAAAPPVARAPEAAGGGSRSEDYSSSPHSAAAAARPLAAEEKQAQSLQPSSSRRSSHYPAAVQSQAAAERGASATAKSRAISILQKKPRHQQLLPSLSSFFFSHRLPDMTAIIKEIVSRNKRRYQEDGFDLDLTYIYPNIIAMGFPAERLEGVYRNNIDDVVRFLDDSKHKNHYKIYNLCAERHYDTAKFNCRVAQYPFEDHNPPQLELIKPFCEDLDQWLSEDDNHVAAIHCKAGKGRTGVMICAYLLHRGKFLKAQEALDFYGEVRTRDKKGVTIPSQRRYVYYYSYLLKNHLDYRPVALLFHKMMFETIPMFSGGTCNPQFVVCQLKVKIYSSNSGPTRREDKFMYFEFPQPLPVCGDIKVEFFHKQNKMLKKDKMFHFWVNTFFIPGPEETSEKVENGSLCDQEIDSICSIERADNDKEYLVLTLTKNDLDKANKDKANRYFSPNFKVKLYFTKTVEEPSNPEASSSTSVTPDVSDNEPDHYRYSDTTDSDPENEPFDEDQHTQITKV*
##    entrezgene
## 1         100
## 2         100
## 3         100
## 4        5728
## 5         100
## 6        5728
## 7        5728
## 8         100
## 9        5728
## 10       5728
## 11       5728

4.10 Retrieve known SNPs located on the human chromosome 8 between positions 148350 and 148612

For this example we’ll first have to connect to a different BioMart database, namely snp.

snpmart = useMart(biomart = "ENSEMBL_MART_SNP", dataset="hsapiens_snp")

The listAttributes() and listFilters() functions give us an overview of the available attributes and filters.
From these we need: refsnp_id, allele, chrom_start and chrom_strand as attributes; and as filters we’ll use: chrom_start, chrom_end and chr_name.
Note that when a chromosome name, a start position and an end position are jointly used as filters, the BioMart webservice interprets this as return everything from the given chromosome between the given start and end positions. Putting our selected attributes and filters into getBM gives:

getBM(attributes = c('refsnp_id','allele','chrom_start','chrom_strand'), 
      filters = c('chr_name','start','end'), 
      values = list(8,148350,148612), 
      mart = snpmart)
##     refsnp_id allele chrom_start chrom_strand
## 1 rs868546642    A/G      148372            1
## 2 rs547420070    A/C      148373            1
## 3  rs77274555    G/A      148391            1
## 4 rs567299969    T/A      148394            1
## 5 rs368076569    G/A      148407            1
## 6 rs745318437    C/G      148497            1
## 7 rs190721891    C/G      148576            1

4.11 Given the human gene TP53, retrieve the human chromosomal location of this gene and also retrieve the chromosomal location and RefSeq id of its homolog in mouse.

The getLDS() (Get Linked Dataset) function provides functionality to link 2 BioMart datasets which each other and construct a query over the two datasets. In Ensembl, linking two datasets translates to retrieving homology data across species. The usage of getLDS is very similar to getBM(). The linked dataset is provided by a separate Mart object and one has to specify filters and attributes for the linked dataset. Filters can either be applied to both datasets or to one of the datasets. Use the listFilters and listAttributes functions on both Mart objects to find the filters and attributes for each dataset (species in Ensembl). The attributes and filters of the linked dataset can be specified with the attributesL and filtersL arguments. Entering all this information into getLDS() gives:

human = useMart("ensembl", dataset = "hsapiens_gene_ensembl")
mouse = useMart("ensembl", dataset = "mmusculus_gene_ensembl")
getLDS(attributes = c("hgnc_symbol","chromosome_name", "start_position"),
       filters = "hgnc_symbol", values = "TP53",mart = human,
      attributesL = c("refseq_mrna","chromosome_name","start_position"), martL = mouse)
##   HGNC.symbol Chromosome.Name Gene.Start..bp. RefSeq.mRNA..e.g..NM_001195597. Chromosome.Name.1 Gene.Start..bp..1
## 1        TP53              17         7661779                                                11          69580359
## 2        TP53              17         7661779                    NM_001127233                11          69580359
## 3        TP53              17         7661779                       NM_011640                11          69580359

5 Using archived versions of Ensembl

It is possible to query archived versions of Ensembl through biomaRt. There are currently two ways to access archived versions.

5.1 Using archive = TRUE

First we list the available Ensembl archives by using the listMarts() function and setting the archive attribute to TRUE. Note that not all archives are available this way and it seems that recently this only gives access to few archives if you don’t see the version of the archive you need please look at the 2nd way to access archives.

listMarts(archive = TRUE)
##                        biomart                     version
## 1              ensembl_mart_51                  Ensembl 51
## 2                  snp_mart_51                      SNP 51
## 3                 vega_mart_51                     Vega 32
## 4              ensembl_mart_50                  Ensembl 50
## 5                  snp_mart_50                      SNP 50
## 6                 vega_mart_50                     Vega 32
## 7              ensembl_mart_49   ENSEMBL GENES 49 (SANGER)
## 8     genomic_features_mart_49            Genomic Features
## 9                  snp_mart_49                         SNP
## 10                vega_mart_49                        Vega
## 11             ensembl_mart_48   ENSEMBL GENES 48 (SANGER)
## 12    genomic_features_mart_48            Genomic Features
## 13                 snp_mart_48                         SNP
## 14                vega_mart_48                        Vega
## 15             ensembl_mart_47   ENSEMBL GENES 47 (SANGER)
## 16    genomic_features_mart_47            Genomic Features
## 17                 snp_mart_47                         SNP
## 18                vega_mart_47                        Vega
## 19    compara_mart_homology_47            Compara homology
## 20 compara_mart_multiple_ga_47 Compara multiple alignments
## 21 compara_mart_pairwise_ga_47 Compara pairwise alignments
## 22             ensembl_mart_46   ENSEMBL GENES 46 (SANGER)
## 23    genomic_features_mart_46            Genomic Features
## 24                 snp_mart_46                         SNP
## 25                vega_mart_46                        Vega
## 26    compara_mart_homology_46            Compara homology
## 27 compara_mart_multiple_ga_46 Compara multiple alignments
## 28 compara_mart_pairwise_ga_46 Compara pairwise alignments
## 29             ensembl_mart_45   ENSEMBL GENES 45 (SANGER)
## 30                 snp_mart_45                         SNP
## 31                vega_mart_45                        Vega
## 32    compara_mart_homology_45            Compara homology
## 33 compara_mart_multiple_ga_45 Compara multiple alignments
## 34 compara_mart_pairwise_ga_45 Compara pairwise alignments
## 35             ensembl_mart_44   ENSEMBL GENES 44 (SANGER)
## 36                 snp_mart_44                         SNP
## 37                vega_mart_44                        Vega
## 38    compara_mart_homology_44            Compara homology
## 39 compara_mart_pairwise_ga_44 Compara pairwise alignments
## 40             ensembl_mart_43   ENSEMBL GENES 43 (SANGER)
## 41                 snp_mart_43                         SNP
## 42                vega_mart_43                        Vega
## 43    compara_mart_homology_43            Compara homology
## 44 compara_mart_pairwise_ga_43 Compara pairwise alignments

Next we select the archive we want to use using the useMart() function, again setting the archive attribute to TRUE and giving the full name of the BioMart e.g. ensembl_mart_46.

ensembl = useMart("ensembl_mart_46", dataset="hsapiens_gene_ensembl", archive = TRUE)
## Note: requested host was redirected from www.ensembl.org to http://aug2007.archive.ensembl.org:80/biomart/martservice
## When using archived Ensembl versions this sometimes can result in connecting to a newer version than the intended Ensembl version
## Check your ensembl version using listMarts(mart)
## Error in useDataset(mart = mart, dataset = dataset, verbose = verbose): The given dataset:  hsapiens_gene_ensembl , is not valid.  Correct dataset names can be obtained with the listDatasets function.

If you don’t know the dataset you want to use could first connect to the BioMart using useMart() and then use the listDatasets() function on this object. After you selected the BioMart database and dataset, queries can be performed in the same way as when using the current BioMart versions.

5.2 Accessing archives through specifying the archive host

Use the http://www.ensembl.org website and go down the bottom of the page. Click on ‘view in Archive’ and select the archive you need. Copy the url and use that url as shown below to connect to the specified BioMart database. The example below shows how to query Ensembl 54.

listMarts(host='may2009.archive.ensembl.org')
##                biomart              version
## 1 ENSEMBL_MART_ENSEMBL           Ensembl 54
## 2     ENSEMBL_MART_SNP Ensembl Variation 54
## 3    ENSEMBL_MART_VEGA              Vega 35
## 4             REACTOME   Reactome(CSHL US) 
## 5     wormbase_current   WormBase (CSHL US)
## 6                pride       PRIDE (EBI UK)
ensembl54=useMart(host='may2009.archive.ensembl.org', biomart='ENSEMBL_MART_ENSEMBL')
ensembl54=useMart(host='may2009.archive.ensembl.org', biomart='ENSEMBL_MART_ENSEMBL', dataset='hsapiens_gene_ensembl')

6 Using a BioMart other than Ensembl

To demonstrate the use of the biomaRt package with non-Ensembl databases the next query is performed using the Wormbase ParaSite BioMart. In this example, we use the listMarts() function to find the name of the available marts, given the URL of Wormbase. We use this to connect to Wormbase BioMart, find and select the gene dataset, and print the first 6 available attributes and filters. Then we use a list of gene names as filter and retrieve associated transcript IDs and the transcript biotype.

listMarts(host = "parasite.wormbase.org")
##         biomart       version
## 1 parasite_mart ParaSite Mart
wormbase = useMart(biomart = "parasite_mart", host = "parasite.wormbase.org")
listDatasets(wormbase)
##     dataset         description version
## 1 wbps_gene All Species (WBPS7)       2
wormbase <- useDataset(mart = wormbase, dataset = "wbps_gene")
head(listFilters(wormbase))
##                  name     description
## 1     species_id_1010          Genome
## 2 nematode_clade_1010  Nematode Clade
## 3     chromosome_name Chromosome name
## 4               start           Start
## 5                 end             End
## 6              strand          Strand
head(listAttributes(wormbase))
##                      name        description         page
## 1          species_id_key      Internal Name feature_page
## 2    production_name_1010     Genome project feature_page
## 3       display_name_1010        Genome name feature_page
## 4        taxonomy_id_1010        Taxonomy ID feature_page
## 5 assembly_accession_1010 Assembly accession feature_page
## 6     nematode_clade_1010     Nematode clade feature_page
getBM(attributes = c("external_gene_id","ensembl_transcript_id","transcript_biotype"), 
      filters="gene_name", 
      values=c("unc-26","his-33"), 
      mart=wormbase)
##   external_gene_id ensembl_transcript_id transcript_biotype
## 1           his-33              F17E9.13     protein_coding
## 2           unc-26               JC8.10a     protein_coding
## 3           unc-26               JC8.10b     protein_coding
## 4           unc-26             JC8.10c.1     protein_coding
## 5           unc-26             JC8.10c.2     protein_coding
## 6           unc-26               JC8.10d     protein_coding

7 biomaRt helper functions

This section describes a set of biomaRt helper functions that can be used to export FASTA format sequences, retrieve values for certain filters and exploring the available filters and attributes in a more systematic manner.

7.1 exportFASTA

The data.frames obtained by the getSequence function can be exported to FASTA files using the exportFASTA() function. One has to specify the data.frame to export and the filename using the file argument.

7.2 Finding out more information on filters

7.2.1 filterType

Boolean filters need a value TRUE or FALSE in biomaRt. Setting the value TRUE will include all information that fulfill the filter requirement. Setting FALSE will exclude the information that fulfills the filter requirement and will return all values that don’t fulfill the filter. For most of the filters, their name indicates if the type is a boolean or not and they will usually start with “with”. However this is not a rule and to make sure you got the type right you can use the function filterType() to investigate the type of the filter you want to use.

filterType("with_affy_hg_u133_plus_2",ensembl)
## [1] "boolean_list"

7.2.2 filterOptions

Some filters have a limited set of values that can be given to them. To know which values these are one can use the filterOptions() function to retrieve the predetermed values of the respective filter.

filterOptions("biotype",ensembl)
## [1] "[3prime_overlapping_ncRNA,antisense,bidirectional_promoter_lncRNA,IG_C_gene,IG_C_pseudogene,IG_D_gene,IG_J_gene,IG_J_pseudogene,IG_pseudogene,IG_V_gene,IG_V_pseudogene,lincRNA,macro_lncRNA,miRNA,misc_RNA,Mt_rRNA,Mt_tRNA,non_coding,polymorphic_pseudogene,processed_pseudogene,processed_transcript,protein_coding,pseudogene,ribozyme,rRNA,scaRNA,scRNA,sense_intronic,sense_overlapping,snoRNA,snRNA,sRNA,TEC,transcribed_processed_pseudogene,transcribed_unitary_pseudogene,transcribed_unprocessed_pseudogene,TR_C_gene,TR_D_gene,TR_J_gene,TR_J_pseudogene,TR_V_gene,TR_V_pseudogene,unitary_pseudogene,unprocessed_pseudogene,vaultRNA]"

If there are no predetermed values e.g. for the entrezgene filter, then filterOptions() will return the type of filter it is. And most of the times the filter name or it’s description will suggest what values one case use for the respective filter (e.g. entrezgene filter will work with enterzgene identifiers as values)

7.3 Attribute Pages

For large BioMart databases such as Ensembl, the number of attributes displayed by the listAttributes() function can be very large. In BioMart databases, attributes are put together in pages, such as sequences, features, homologs for Ensembl. An overview of the attributes pages present in the respective BioMart dataset can be obtained with the attributePages() function.

pages = attributePages(ensembl)
pages
## [1] "feature_page" "structure"    "homologs"     "snp"          "snp_somatic"  "sequences"

To show us a smaller list of attributes which belong to a specific page, we can now specify this in the listAttributes() function. The set of attributes is still quite long, so we use head() to show only the first few items here.

head(listAttributes(ensembl, page="feature_page"))
##                    name           description         page
## 1       ensembl_gene_id       Ensembl Gene ID feature_page
## 2 ensembl_transcript_id Ensembl Transcript ID feature_page
## 3    ensembl_peptide_id    Ensembl Protein ID feature_page
## 4       ensembl_exon_id       Ensembl Exon ID feature_page
## 5           description           Description feature_page
## 6       chromosome_name       Chromosome Name feature_page

We now get a short list of attributes related to the region where the genes are located.

8 Local BioMart databases

The biomaRt package can be used with a local install of a public BioMart database or a locally developed BioMart database and web service. In order for biomaRt to recognize the database as a BioMart, make sure that the local database you create has a name conform with database_mart_version where database is the name of the database and version is a version number. No more underscores than the ones showed should be present in this name. A possible name is for example ensemblLocal_mart_46. ## Minimum requirements for local database installation More information on installing a local copy of a BioMart database or develop your own BioMart database and webservice can be found on http://www.biomart.org Once the local database is installed you can use biomaRt on this database by:

listMarts(host="www.myLocalHost.org", path="/myPathToWebservice/martservice")
mart=useMart("nameOfMyMart",dataset="nameOfMyDataset",host="www.myLocalHost.org", path="/myPathToWebservice/martservice")

For more information on how to install a public BioMart database see: http://www.biomart.org/install.html and follow link databases.

9 Using select()

In order to provide a more consistent interface to all annotations in Bioconductor the select(), columns(), keytypes() and keys() have been implemented to wrap some of the existing functionality above. These methods can be called in the same manner that they are used in other parts of the project except that instead of taking a AnnotationDb derived class they take instead a Mart derived class as their 1st argument. Otherwise usage should be essentially the same. You still use columns() to discover things that can be extracted from a Mart, and keytypes() to discover which things can be used as keys with select().

mart <- useMart(dataset="hsapiens_gene_ensembl",biomart='ensembl')
head(keytypes(mart), n=3)
## [1] "affy_hc_g110"        "affy_hg_focus"       "affy_hg_u133_plus_2"
head(columns(mart), n=3)
## [1] "3_utr_end"   "3_utr_end"   "3_utr_start"

And you still can use keys() to extract potential keys, for a particular key type.

k = keys(mart, keytype="chromosome_name")
head(k, n=3)
## [1] "1" "2" "3"

When using keys(), you can even take advantage of the extra arguments that are available for others keys methods.

k = keys(mart, keytype="chromosome_name", pattern="LRG")
head(k, n=3)
## character(0)

Unfortunately the keys() method will not work with all key types because they are not all supported.

But you can still use select() here to extract columns of data that match a particular set of keys (this is basically a wrapper for getBM()).

affy=c("202763_at","209310_s_at","207500_at")
select(mart, keys=affy, columns=c('affy_hg_u133_plus_2','entrezgene'),
  keytype='affy_hg_u133_plus_2')
##   affy_hg_u133_plus_2 entrezgene
## 1           207500_at        838
## 2           202763_at        836
## 3         209310_s_at        837

So why would we want to do this when we already have functions like getBM()? For two reasons: 1) for people who are familiar with select and it’s helper methods, they can now proceed to use biomaRt making the same kinds of calls that are already familiar to them and 2) because the select method is implemented in many places elsewhere, the fact that these methods are shared allows for more convenient programmatic access of all these resources. An example of a package that takes advantage of this is the OrganismDbi package. Where several packages can be accessed as if they were one resource.

10 Session Info

sessionInfo()
## R version 3.3.1 (2016-06-21)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 16.04.1 LTS
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C               LC_TIME=en_US.UTF-8        LC_COLLATE=C              
##  [5] LC_MONETARY=en_US.UTF-8    LC_MESSAGES=en_US.UTF-8    LC_PAPER=en_US.UTF-8       LC_NAME=C                 
##  [9] LC_ADDRESS=C               LC_TELEPHONE=C             LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C       
## 
## attached base packages:
## [1] stats     graphics  grDevices utils     datasets  methods   base     
## 
## other attached packages:
## [1] biomaRt_2.30.0  BiocStyle_2.2.0
## 
## loaded via a namespace (and not attached):
##  [1] Rcpp_0.12.7          IRanges_2.8.0        XML_3.98-1.4         digest_0.6.10        assertthat_0.1      
##  [6] bitops_1.0-6         DBI_0.5-1            stats4_3.3.1         formatR_1.4          magrittr_1.5        
## [11] RSQLite_1.0.0        evaluate_0.10        stringi_1.1.2        S4Vectors_0.12.0     rmarkdown_1.1       
## [16] tools_3.3.1          stringr_1.1.0        Biobase_2.34.0       RCurl_1.95-4.8       parallel_3.3.1      
## [21] yaml_2.1.13          BiocGenerics_0.20.0  AnnotationDbi_1.36.0 htmltools_0.3.5      knitr_1.14          
## [26] tibble_1.2
warnings()
## NULL