Ensembl-REST

https://img.shields.io/badge/Say%20Thanks-!-1EAEDB.svg

A Python interface to the Ensembl REST APIs. A whole world of biological data at your fingertips.

The Ensembl database contains reference biological data on almost any organism. Now it is easy to access this data programatically through their REST API.

The full list of endpoints for the Ensembl REST API endpoints along with endpoint-specific documentation can be found on their website.

This library also includes some utilities built on top of the APIs designed to ease working with them, including an AssemblyMapper class that helps in the conversion between different genome assemblies.

This project uses code from RESTEasy, which made my life much easier. Thanks!

Installation

You can install from PyPI:

$ pip install ensembl_rest

Examples

The library exports two main classes: ensembl_rest.EnsemblClient and ensembl_rest.EnsemblGenomesClient that point respectively to the main REST API and to the Ensembl Genomes REST API.

>>> import ensembl_rest

>>> client = ensembl_rest.EnsemblClient()

If you want to use a method from the REST API, say: GET lookup/symbol/:species/:symbol (http://rest.ensembl.org/documentation/info/symbol_lookup) then the corresponding method on the class is called after the last string in the link to the documentation page, in this case, symbol_lookup.

>>> help(client.symbol_lookup)
Help on method symbol_lookup in module ensembl_rest.core.baseclient:

symbol_lookup(*args, **kwargs) method of ensembl_rest.core.baseclient.EnsemblClient instance
    ``GET lookup/symbol/:species/:symbol``

    Find the species and database for a symbol in a linked external database
    - More info: https://rest.ensembl.org/documentation/info/symbol_lookup

We can see from the resource string GET lookup/symbol/:species/:symbol that this method contains 2 parameters called species and symbol, so we can call the method in the following way:

>>> client.symbol_lookup(species='homo sapiens',
                         symbol='TP53')

# Or like this...
>>> client.symbol_lookup('homo sapiens', 'TP53')
{'source': 'ensembl_havana',
  'object_type': 'Gene',
  'logic_name': 'ensembl_havana_gene',
  'version': 14,
  'species': 'human',
  'description': 'BRCA2, DNA repair associated [Source:HGNC Symbol;Acc:HGNC:1101]',
  'display_name': 'BRCA2',
  'assembly_name': 'GRCh38',
  'biotype': 'protein_coding',
  'end': 32400266,
  'seq_region_name': '13',
  'db_type': 'core',
  'strand': 1,
  'id': 'ENSG00000139618',
  'start': 32315474}

One can pass request parameters to modify the response with the params keyword (the specific parameters to pass depend on the specific endpoint, the official endpoints documentation can be found here, the parameters to be passed this way are the optional ones for each endpoint):

>>> client.symbol_lookup('human', 'BRCA2', params={'expand':True})
{'source': 'ensembl_havana',
 'object_type': 'Gene',
 'logic_name': 'ensembl_havana_gene',
 'seq_region_name': '13',
 'db_type': 'core',
 'strand': 1,
 'id': 'ENSG00000139618',
 'Transcript': [{'source': 'ensembl_havana',
   'object_type': 'Transcript',
   'logic_name': 'ensembl_havana_transcript',
   'Exon': [{'object_type': 'Exon',
     'version': 4,
     'species': 'human',
     'assembly_name': 'GRCh38',
     ...
     ...
     ...
 'biotype': 'protein_coding',
 'start': 32315474}

In that way, you can pass the parameters for the POST endpoints, such as:

>>> client.symbol_post(species='human',
                       params={'symbols': ["BRCA2",
                                           "TP53",
                                           "BRAF" ]})
{
    "BRCA2": {
        "source": "ensembl_havana",
        "object_type": "Gene",
        "logic_name": "ensembl_havana_gene",
        "version": 14,
        "species": "homo_sapiens",
        "description": "BRCA2, DNA repair associated [Source:HGNC Symbol;Acc:HGNC:1101]",
        ...
    },
    "TP53": {
        ...
    }.
    "BRAF": {
        ...
        "strand": -1,
        "id": "ENSG00000157764",
        "start": 140719327
    }
}

Another common usage is to fetch sequences of known genes:

>>> client.sequence_id('ENSG00000157764')
{'desc': 'chromosome:GRCh38:7:140719327:140924928:-1',
 'query': 'ENSG00000157764',
 'version': 13,
 'id': 'ENSG00000157764',
 'seq': 'TTCCCCCAATCCCCTCAGGCTCGGCTGCGCCCGGGGC...ACTGCTATAATAAAGATTGACTGCATGGAGAAGTCTTCA',
 'molecule': 'dna'}

Or to map betweeen assemblies…

>>> client.assembly_map(species='human',
                        asm_one='GRCh37',
                        region='X:1000000..1000100:1',
                        asm_two='GRCh38')


# Or...
>>> region_str = ensembl_rest.region_str(chom='X',
                                         start=1000000,
                                         end=1000100)

>>> client.assembly_map(species='human',
                        asm_one='GRCh37',
                        region=region_str,
                        asm_two='GRCh38')
{'mappings': [{'original': {'seq_region_name': 'X',
    'strand': 1,
    'coord_system': 'chromosome',
    'end': 1000100,
    'start': 1000000,
    'assembly': 'GRCh37'},
   'mapped': {'seq_region_name': 'X',
    'strand': 1,
    'coord_system': 'chromosome',
    'end': 1039365,
    'start': 1039265,
    'assembly': 'GRCh38'}}]}

The above problem (mapping from one assembly to another) is so frequent that the library provides a specialized class AssemblyMapper to efficiently mapping large amounts of regions between assemblies. This class avoids the time-consuming task of making a web request every time a mapping is needed by fetching the mapping of the whole assembly right from the instantiation. This is a time-consuming operation by itself, but it pays off when one has to transform repeatedly betweeen assemblies.:

>>> mapper = ensembl_rest.AssemblyMapper(from_assembly='GRCh37'
...                                      to_assembly='GRCh38')

>>> mapper.map(chrom='1', pos=1000000)
1064620

Meta

Author: Ad115 - Githuba.garcia230395@gmail.com

Project pages: Docs - @GitHub - @PyPI

Distributed under the MIT license. See LICENSE for more information.

Contributing

  1. Check for open issues or open a fresh issue to start a discussion around a feature idea or a bug.
  2. Fork the repository on GitHub to start making your changes to a feature branch, derived from the master branch.
  3. Write a test which shows that the bug was fixed or that the feature works as expected.
  4. Send a pull request and bug the maintainer until it gets merged and published.

API Documentation

ensembl_rest.region_str(chrom, start, end=None, strand=1)[source]

Assemble a region string suitable for consumption for the Ensembl REST API.

The generated string has the format: {chrom}:{start}..{end}:{strand}

class ensembl_rest.AssemblyMapper(from_assembly='GRCh37', to_assembly='GRCh38', species='human')[source]

Structure that optimizes the mapping between diferent genome assemblies.

Example:

>>> mapper = AssemblyMapper(from_assembly='GRCh37'
...                         to_assembly='GRCh38')

>>> mapper.map(chrom='1', pos=1000000)
1064620
map(chrom, pos)[source]

Map the given position.

The mapping is between the specified assemblies when creating the object. (default: map position from assembly GRCh37 to assembly GRCh38)

class ensembl_rest.EnsemblClient(base_url=None)

A client for the Ensembl REST API (https://rest.ensembl.org/)

VariantAnnotationSet(*args, **kwargs)

POST ga4gh/variantannotationsets/search

Return a list of annotation sets in GA4GH format - More info: https://rest.ensembl.org/documentation/info/VariantAnnotationSet

VariantAnnotationSet_id(*args, **kwargs)

GET ga4gh/variantannotationsets/:id

Return meta data for a specific annotation set in GA4GH format - More info: https://rest.ensembl.org/documentation/info/VariantAnnotationSet_id

analysis(*args, **kwargs)

GET info/analysis/:species

List the names of analyses involved in generating Ensembl data. - More info: https://rest.ensembl.org/documentation/info/analysis

archive_id_get(*args, **kwargs)

GET archive/id/:id

Uses the given identifier to return the archived sequence - More info: https://rest.ensembl.org/documentation/info/archive_id_get

archive_id_post(*args, **kwargs)

POST archive/id

Retrieve the archived sequence for a set of identifiers - More info: https://rest.ensembl.org/documentation/info/archive_id_post

array(*args, **kwargs)

GET regulatory/species/:species/microarray/:microarray/vendor/:vendor

Returns information about a specific microarray - More info: https://rest.ensembl.org/documentation/info/array

assembly_cdna(*args, **kwargs)

GET map/cdna/:id/:region

Convert from cDNA coordinates to genomic coordinates. Output reflects forward orientation coordinates as returned from the Ensembl API. - More info: https://rest.ensembl.org/documentation/info/assembly_cdna

assembly_cds(*args, **kwargs)

GET map/cds/:id/:region

Convert from CDS coordinates to genomic coordinates. Output reflects forward orientation coordinates as returned from the Ensembl API. - More info: https://rest.ensembl.org/documentation/info/assembly_cds

assembly_info(*args, **kwargs)

GET info/assembly/:species

List the currently available assemblies for a species, along with toplevel sequences, chromosomes and cytogenetic bands. - More info: https://rest.ensembl.org/documentation/info/assembly_info

assembly_map(*args, **kwargs)

GET map/:species/:asm_one/:region/:asm_two

Convert the co-ordinates of one assembly to another - More info: https://rest.ensembl.org/documentation/info/assembly_map

assembly_stats(*args, **kwargs)

GET info/assembly/:species/:region_name

Returns information about the specified toplevel sequence region for the given species. - More info: https://rest.ensembl.org/documentation/info/assembly_stats

assembly_translation(*args, **kwargs)

GET map/translation/:id/:region

Convert from protein (translation) coordinates to genomic coordinates. Output reflects forward orientation coordinates as returned from the Ensembl API. - More info: https://rest.ensembl.org/documentation/info/assembly_translation

beacon_get(*args, **kwargs)

GET ga4gh/beacon

Return Beacon information - More info: https://rest.ensembl.org/documentation/info/beacon_get

beacon_query_get(*args, **kwargs)

GET ga4gh/beacon/query

Return the Beacon response for allele information - More info: https://rest.ensembl.org/documentation/info/beacon_query_get

beacon_query_post(*args, **kwargs)

POST ga4gh/beacon/query

Return the Beacon response for allele information - More info: https://rest.ensembl.org/documentation/info/beacon_query_post

biotypes(*args, **kwargs)

GET info/biotypes/:species

List the functional classifications of gene models that Ensembl associates with a particular species. Useful for restricting the type of genes/transcripts retrieved by other endpoints. - More info: https://rest.ensembl.org/documentation/info/biotypes

biotypes_groups(*args, **kwargs)

GET info/biotypes/groups/:group/:object_type

Without argument the list of available biotype groups is returned. With :group argument provided, list the properties of biotypes within that group. Object type (gene or transcript) can be provided for filtering. - More info: https://rest.ensembl.org/documentation/info/biotypes_groups

biotypes_name(*args, **kwargs)

GET info/biotypes/name/:name/:object_type

List the properties of biotypes with a given name. Object type (gene or transcript) can be provided for filtering. - More info: https://rest.ensembl.org/documentation/info/biotypes_name

cafe_tree(*args, **kwargs)

GET cafe/genetree/id/:id

Retrieves a cafe tree of the gene tree using the gene tree stable identifier - More info: https://rest.ensembl.org/documentation/info/cafe_tree

cafe_tree_member_id(*args, **kwargs)

GET cafe/genetree/member/id/:id

Retrieves the cafe tree of the gene tree that contains the gene / transcript / translation stable identifier - More info: https://rest.ensembl.org/documentation/info/cafe_tree_member_id

cafe_tree_member_symbol(*args, **kwargs)

GET cafe/genetree/member/symbol/:species/:symbol

Retrieves the cafe tree of the gene tree that contains the gene identified by a symbol - More info: https://rest.ensembl.org/documentation/info/cafe_tree_member_symbol

compara_methods(*args, **kwargs)

GET info/compara/methods

List all compara analyses available (an analysis defines the type of comparative data). - More info: https://rest.ensembl.org/documentation/info/compara_methods

compara_species_sets(*args, **kwargs)

GET info/compara/species_sets/:method

List all collections of species analysed with the specified compara method. - More info: https://rest.ensembl.org/documentation/info/compara_species_sets

comparas(*args, **kwargs)

GET info/comparas

Lists all available comparative genomics databases and their data release. - More info: https://rest.ensembl.org/documentation/info/comparas

data(*args, **kwargs)

GET info/data

Shows the data releases available on this REST server. May return more than one release (unfrequent non-standard Ensembl configuration). - More info: https://rest.ensembl.org/documentation/info/data

external_dbs(*args, **kwargs)

GET info/external_dbs/:species

Lists all available external sources for a species. - More info: https://rest.ensembl.org/documentation/info/external_dbs

family(*args, **kwargs)

GET family/id/:id

Retrieves a family information using the family stable identifier - More info: https://rest.ensembl.org/documentation/info/family

family_member_id(*args, **kwargs)

GET family/member/id/:id

Retrieves the information for all the families that contains the gene / transcript / translation stable identifier - More info: https://rest.ensembl.org/documentation/info/family_member_id

family_member_symbol(*args, **kwargs)

GET family/member/symbol/:species/:symbol

Retrieves the information for all the families that contains the gene identified by a symbol - More info: https://rest.ensembl.org/documentation/info/family_member_symbol

features_id(*args, **kwargs)

GET ga4gh/features/:id

Return the GA4GH record for a specific sequence feature given its identifier - More info: https://rest.ensembl.org/documentation/info/features_id

features_post(*args, **kwargs)

POST ga4gh/features/search

Return a list of sequence annotation features in GA4GH format - More info: https://rest.ensembl.org/documentation/info/features_post

fetch_all_epigenomes(*args, **kwargs)

GET regulatory/species/:species/epigenome

Returns information about all epigenomes available for the given species - More info: https://rest.ensembl.org/documentation/info/fetch_all_epigenomes

gacallSet(*args, **kwargs)

POST ga4gh/callsets/search

Return a list of sets of genotype calls for specific samples in GA4GH format - More info: https://rest.ensembl.org/documentation/info/gacallSet

gacallset_id(*args, **kwargs)

GET ga4gh/callsets/:id

Return the GA4GH record for a specific CallSet given its identifier - More info: https://rest.ensembl.org/documentation/info/gacallset_id

gadataset(*args, **kwargs)

POST ga4gh/datasets/search

Return a list of datasets in GA4GH format - More info: https://rest.ensembl.org/documentation/info/gadataset

gadataset_id(*args, **kwargs)

GET ga4gh/datasets/:id

Return the GA4GH record for a specific dataset given its identifier - More info: https://rest.ensembl.org/documentation/info/gadataset_id

gafeatureset(*args, **kwargs)

POST ga4gh/featuresets/search

Return a list of feature sets in GA4GH format - More info: https://rest.ensembl.org/documentation/info/gafeatureset

gafeatureset_id(*args, **kwargs)

GET ga4gh/featuresets/:id

Return the GA4GH record for a specific featureSet given its identifier - More info: https://rest.ensembl.org/documentation/info/gafeatureset_id

gavariant_id(*args, **kwargs)

GET ga4gh/variants/:id

Return the GA4GH record for a specific variant given its identifier. - More info: https://rest.ensembl.org/documentation/info/gavariant_id

gavariantannotations(*args, **kwargs)

POST ga4gh/variantannotations/search

Return variant annotation information in GA4GH format for a region on a reference sequence - More info: https://rest.ensembl.org/documentation/info/gavariantannotations

gavariants(*args, **kwargs)

POST ga4gh/variants/search

Return variant call information in GA4GH format for a region on a reference sequence - More info: https://rest.ensembl.org/documentation/info/gavariants

gavariantset(*args, **kwargs)

POST ga4gh/variantsets/search

Return a list of variant sets in GA4GH format - More info: https://rest.ensembl.org/documentation/info/gavariantset

gavariantset_id(*args, **kwargs)

GET ga4gh/variantsets/:id

Return the GA4GH record for a specific VariantSet given its identifier - More info: https://rest.ensembl.org/documentation/info/gavariantset_id

genetree(*args, **kwargs)

GET genetree/id/:id

Retrieves a gene tree for a gene tree stable identifier - More info: https://rest.ensembl.org/documentation/info/genetree

genetree_member_id(*args, **kwargs)

GET genetree/member/id/:id

Retrieves the gene tree that contains the gene / transcript / translation stable identifier - More info: https://rest.ensembl.org/documentation/info/genetree_member_id

genetree_member_symbol(*args, **kwargs)

GET genetree/member/symbol/:species/:symbol

Retrieves the gene tree that contains the gene identified by a symbol - More info: https://rest.ensembl.org/documentation/info/genetree_member_symbol

genomic_alignment_region(*args, **kwargs)

GET alignment/region/:species/:region

Retrieves genomic alignments as separate blocks based on a region and species - More info: https://rest.ensembl.org/documentation/info/genomic_alignment_region

homology_ensemblgene(*args, **kwargs)

GET homology/id/:id

Retrieves homology information (orthologs) by Ensembl gene id - More info: https://rest.ensembl.org/documentation/info/homology_ensemblgene

homology_symbol(*args, **kwargs)

GET homology/symbol/:species/:symbol

Retrieves homology information (orthologs) by symbol - More info: https://rest.ensembl.org/documentation/info/homology_symbol

ld_id_get(*args, **kwargs)

GET ld/:species/:id/:population_name

Computes and returns LD values between the given variant and all other variants in a window centered around the given variant. The window size is set to 500 kb. - More info: https://rest.ensembl.org/documentation/info/ld_id_get

ld_pairwise_get(*args, **kwargs)

GET ld/:species/pairwise/:id1/:id2

Computes and returns LD values between the given variants. - More info: https://rest.ensembl.org/documentation/info/ld_pairwise_get

ld_region_get(*args, **kwargs)

GET ld/:species/region/:region/:population_name

Computes and returns LD values between all pairs of variants in the defined region. - More info: https://rest.ensembl.org/documentation/info/ld_region_get

list_all_microarrays(*args, **kwargs)

GET regulatory/species/:species/microarray

Returns information about all microarrays available for the given species - More info: https://rest.ensembl.org/documentation/info/list_all_microarrays

lookup(*args, **kwargs)

GET lookup/id/:id

Find the species and database for a single identifier e.g. gene, transcript, protein - More info: https://rest.ensembl.org/documentation/info/lookup

lookup_post(*args, **kwargs)

POST lookup/id

Find the species and database for several identifiers. IDs that are not found are returned with no data. - More info: https://rest.ensembl.org/documentation/info/lookup_post

ontology_ancestors(*args, **kwargs)

GET ontology/ancestors/:id

Reconstruct the entire ancestry of a term from is_a and part_of relationships - More info: https://rest.ensembl.org/documentation/info/ontology_ancestors

ontology_ancestors_chart(*args, **kwargs)

GET ontology/ancestors/chart/:id

Reconstruct the entire ancestry of a term from is_a and part_of relationships. - More info: https://rest.ensembl.org/documentation/info/ontology_ancestors_chart

ontology_descendants(*args, **kwargs)

GET ontology/descendants/:id

Find all the terms descended from a given term. By default searches are conducted within the namespace of the given identifier - More info: https://rest.ensembl.org/documentation/info/ontology_descendants

ontology_id(*args, **kwargs)

GET ontology/id/:id

Search for an ontological term by its namespaced identifier - More info: https://rest.ensembl.org/documentation/info/ontology_id

ontology_name(*args, **kwargs)

GET ontology/name/:name

Search for a list of ontological terms by their name - More info: https://rest.ensembl.org/documentation/info/ontology_name

overlap_id(*args, **kwargs)

GET overlap/id/:id

Retrieves features (e.g. genes, transcripts, variants and more) that overlap a region defined by the given identifier. - More info: https://rest.ensembl.org/documentation/info/overlap_id

overlap_region(*args, **kwargs)

GET overlap/region/:species/:region

Retrieves features (e.g. genes, transcripts, variants and more) that overlap a given region. - More info: https://rest.ensembl.org/documentation/info/overlap_region

overlap_translation(*args, **kwargs)

GET overlap/translation/:id

Retrieve features related to a specific Translation as described by its stable ID (e.g. domains, variants). - More info: https://rest.ensembl.org/documentation/info/overlap_translation

phenotype_accession(*args, **kwargs)

GET /phenotype/accession/:species/:accession

Return phenotype annotations for genomic features given a phenotype ontology accession - More info: https://rest.ensembl.org/documentation/info/phenotype_accession

phenotype_gene(*args, **kwargs)

GET /phenotype/gene/:species/:gene

Return phenotype annotations for a given gene. - More info: https://rest.ensembl.org/documentation/info/phenotype_gene

phenotype_region(*args, **kwargs)

GET /phenotype/region/:species/:region

Return phenotype annotations that overlap a given genomic region. - More info: https://rest.ensembl.org/documentation/info/phenotype_region

phenotype_term(*args, **kwargs)

GET /phenotype/term/:species/:term

Return phenotype annotations for genomic features given a phenotype ontology term - More info: https://rest.ensembl.org/documentation/info/phenotype_term

ping(*args, **kwargs)

GET info/ping

Checks if the service is alive. - More info: https://rest.ensembl.org/documentation/info/ping

probe(*args, **kwargs)

GET regulatory/species/:species/microarray/:microarray/probe/:probe

Returns information about a specific probe from a microarray - More info: https://rest.ensembl.org/documentation/info/probe

probe_set(*args, **kwargs)

GET regulatory/species/:species/microarray/:microarray/probe_set/:probe_set

Returns information about a specific probe_set from a microarray - More info: https://rest.ensembl.org/documentation/info/probe_set

referenceSets(*args, **kwargs)

POST ga4gh/referencesets/search

Return a list of reference sets in GA4GH format - More info: https://rest.ensembl.org/documentation/info/referenceSets

referenceSets_id(*args, **kwargs)

GET ga4gh/referencesets/:id

Return data for a specific reference set in GA4GH format - More info: https://rest.ensembl.org/documentation/info/referenceSets_id

references(*args, **kwargs)

POST ga4gh/references/search

Return a list of reference sequences in GA4GH format - More info: https://rest.ensembl.org/documentation/info/references

references_id(*args, **kwargs)

GET ga4gh/references/:id

Return data for a specific reference in GA4GH format by id - More info: https://rest.ensembl.org/documentation/info/references_id

regulatory_id(*args, **kwargs)

GET regulatory/species/:species/id/:id

Returns a RegulatoryFeature given its stable ID (e.g. ENSR00000099113) - More info: https://rest.ensembl.org/documentation/info/regulatory_id

rest(*args, **kwargs)

GET info/rest

Shows the current version of the Ensembl REST API. - More info: https://rest.ensembl.org/documentation/info/rest

sequence_id(*args, **kwargs)

GET sequence/id/:id

Request multiple types of sequence by stable identifier. Supports feature masking and expand options. - More info: https://rest.ensembl.org/documentation/info/sequence_id

sequence_id_post(*args, **kwargs)

POST sequence/id

Request multiple types of sequence by a stable identifier list. - More info: https://rest.ensembl.org/documentation/info/sequence_id_post

sequence_region(*args, **kwargs)

GET sequence/region/:species/:region

Returns the genomic sequence of the specified region of the given species. Supports feature masking and expand options. - More info: https://rest.ensembl.org/documentation/info/sequence_region

sequence_region_post(*args, **kwargs)

POST sequence/region/:species

Request multiple types of sequence by a list of regions. - More info: https://rest.ensembl.org/documentation/info/sequence_region_post

software(*args, **kwargs)

GET info/software

Shows the current version of the Ensembl API used by the REST server. - More info: https://rest.ensembl.org/documentation/info/software

species(*args, **kwargs)

GET info/species

Lists all available species, their aliases, available adaptor groups and data release. - More info: https://rest.ensembl.org/documentation/info/species

species_id(*args, **kwargs)

GET eqtl/stable_id/:species/:stable_id

Returns the p-value for each SNP in a given gene (e.g. ENSG00000227232) - More info: https://rest.ensembl.org/documentation/info/species_id

species_variant(*args, **kwargs)

GET eqtl/variant_name/:species/:variant_name

Returns the p-values for a SNP (e.g. rs123) - More info: https://rest.ensembl.org/documentation/info/species_variant

symbol_lookup(*args, **kwargs)

GET lookup/symbol/:species/:symbol

Find the species and database for a symbol in a linked external database - More info: https://rest.ensembl.org/documentation/info/symbol_lookup

symbol_post(*args, **kwargs)

POST lookup/symbol/:species/:symbol

Find the species and database for a set of symbols in a linked external database. Unknown symbols are omitted from the response. - More info: https://rest.ensembl.org/documentation/info/symbol_post

taxonomy_classification(*args, **kwargs)

GET taxonomy/classification/:id

Return the taxonomic classification of a taxon node - More info: https://rest.ensembl.org/documentation/info/taxonomy_classification

taxonomy_id(*args, **kwargs)

GET taxonomy/id/:id

Search for a taxonomic term by its identifier or name - More info: https://rest.ensembl.org/documentation/info/taxonomy_id

taxonomy_name(*args, **kwargs)

GET taxonomy/name/:name

Search for a taxonomic id by a non-scientific name - More info: https://rest.ensembl.org/documentation/info/taxonomy_name

tissues(*args, **kwargs)

GET eqtl/tissue/:species/

Returns all tissues currently available in the DB - More info: https://rest.ensembl.org/documentation/info/tissues

transcript_haplotypes_get(*args, **kwargs)

GET transcript_haplotypes/:species/:id

Computes observed transcript haplotype sequences based on phased genotype data - More info: https://rest.ensembl.org/documentation/info/transcript_haplotypes_get

variant_recoder(*args, **kwargs)

GET variant_recoder/:species/:id

Translate a variant identifier or HGVS notation to all possible variant IDs and HGVS - More info: https://rest.ensembl.org/documentation/info/variant_recoder

variant_recoder_post(*args, **kwargs)

POST variant_recoder/:species

Translate a list of variant identifiers or HGVS notations to all possible variant IDs and HGVS - More info: https://rest.ensembl.org/documentation/info/variant_recoder_post

variation(*args, **kwargs)

GET info/variation/:species

List the variation sources used in Ensembl for a species. - More info: https://rest.ensembl.org/documentation/info/variation

variation_consequence_types(*args, **kwargs)

GET info/variation/consequence_types

Lists all variant consequence types. - More info: https://rest.ensembl.org/documentation/info/variation_consequence_types

variation_id(*args, **kwargs)

GET variation/:species/:id

Uses a variant identifier (e.g. rsID) to return the variation features including optional genotype, phenotype and population data - More info: https://rest.ensembl.org/documentation/info/variation_id

variation_pmcid_get(*args, **kwargs)

GET variation/:species/pmcid/:pmcid

Fetch variants by publication using PubMed Central reference number (PMCID) - More info: https://rest.ensembl.org/documentation/info/variation_pmcid_get

variation_pmid_get(*args, **kwargs)

GET variation/:species/pmid/:pmid

Fetch variants by publication using PubMed reference number (PMID) - More info: https://rest.ensembl.org/documentation/info/variation_pmid_get

variation_populations(*args, **kwargs)

GET info/variation/populations/:species

List all populations for a species - More info: https://rest.ensembl.org/documentation/info/variation_populations

variation_post(*args, **kwargs)

POST variation/:species/

Uses a list of variant identifiers (e.g. rsID) to return the variation features including optional genotype, phenotype and population data - More info: https://rest.ensembl.org/documentation/info/variation_post

vep_hgvs_get(*args, **kwargs)

GET vep/:species/hgvs/:hgvs_notation

Fetch variant consequences based on a HGVS notation - More info: https://rest.ensembl.org/documentation/info/vep_hgvs_get

vep_hgvs_post(*args, **kwargs)

POST vep/:species/hgvs

Fetch variant consequences for multiple HGVS notations - More info: https://rest.ensembl.org/documentation/info/vep_hgvs_post

vep_id_get(*args, **kwargs)

GET vep/:species/id/:id

Fetch variant consequences based on a variant identifier - More info: https://rest.ensembl.org/documentation/info/vep_id_get

vep_id_post(*args, **kwargs)

POST vep/:species/id

Fetch variant consequences for multiple ids - More info: https://rest.ensembl.org/documentation/info/vep_id_post

vep_region_get(*args, **kwargs)

GET vep/:species/region/:region/:allele/

Fetch variant consequences - More info: https://rest.ensembl.org/documentation/info/vep_region_get

vep_region_post(*args, **kwargs)

POST vep/:species/region

Fetch variant consequences for multiple regions - More info: https://rest.ensembl.org/documentation/info/vep_region_post

xref_external(*args, **kwargs)

GET xrefs/symbol/:species/:symbol

Looks up an external symbol and returns all Ensembl objects linked to it. This can be a display name for a gene/transcript/translation, a synonym or an externally linked reference. If a gene’s transcript is linked to the supplied symbol the service will return both gene and transcript (it supports transient links). - More info: https://rest.ensembl.org/documentation/info/xref_external

xref_id(*args, **kwargs)

GET xrefs/id/:id

Perform lookups of Ensembl Identifiers and retrieve their external references in other databases - More info: https://rest.ensembl.org/documentation/info/xref_id

xref_name(*args, **kwargs)

GET xrefs/name/:species/:name

Performs a lookup based upon the primary accession or display label of an external reference and returning the information we hold about the entry - More info: https://rest.ensembl.org/documentation/info/xref_name

class ensembl_rest.EnsemblGenomesClient(base_url=None)

A client for the EnsemblGenomes REST API (http://rest.ensemblgenomes.org)

analysis(*args, **kwargs)

GET info/analysis/:species

List the names of analyses involved in generating Ensembl data. - More info: /documentation/info/analysis

array_info(*args, **kwargs)

GET regulatory/species/:species/microarray/:microarray/vendor/:vendor

Returns information about a specific microarray - More info: /documentation/info/array_info

assembly_cdna(*args, **kwargs)

GET map/cdna/:id/:region

Convert from cDNA coordinates to genomic coordinates. Output reflects forward orientation coordinates as returned from the Ensembl API. - More info: /documentation/info/assembly_cdna

assembly_cds(*args, **kwargs)

GET map/cds/:id/:region

Convert from CDS coordinates to genomic coordinates. Output reflects forward orientation coordinates as returned from the Ensembl API. - More info: /documentation/info/assembly_cds

assembly_info(*args, **kwargs)

GET info/assembly/:species

List the currently available assemblies for a species, along with toplevel sequences, chromosomes and cytogenetic bands. - More info: /documentation/info/assembly_info

assembly_map(*args, **kwargs)

GET map/:species/:asm_one/:region/:asm_two

Convert the co-ordinates of one assembly to another - More info: /documentation/info/assembly_map

assembly_stats(*args, **kwargs)

GET info/assembly/:species/:region_name

Returns information about the specified toplevel sequence region for the given species. - More info: /documentation/info/assembly_stats

assembly_translation(*args, **kwargs)

GET map/translation/:id/:region

Convert from protein (translation) coordinates to genomic coordinates. Output reflects forward orientation coordinates as returned from the Ensembl API. - More info: /documentation/info/assembly_translation

biotypes(*args, **kwargs)

GET info/biotypes/:species

List the functional classifications of gene models that Ensembl associates with a particular species. Useful for restricting the type of genes/transcripts retrieved by other endpoints. - More info: /documentation/info/biotypes

biotypes_groups(*args, **kwargs)

GET info/biotypes/groups/:group/:object_type

Without argument the list of available biotype groups is returned. With :group argument provided, list the properties of biotypes within that group. Object type (gene or transcript) can be provided for filtering. - More info: /documentation/info/biotypes_groups

biotypes_name(*args, **kwargs)

GET info/biotypes/name/:name/:object_type

List the properties of biotypes with a given name. Object type (gene or transcript) can be provided for filtering. - More info: /documentation/info/biotypes_name

compara_methods(*args, **kwargs)

GET info/compara/methods

List all compara analyses available (an analysis defines the type of comparative data). - More info: /documentation/info/compara_methods

compara_species_sets(*args, **kwargs)

GET info/compara/species_sets/:method

List all collections of species analysed with the specified compara method. - More info: /documentation/info/compara_species_sets

comparas(*args, **kwargs)

GET info/comparas

Lists all available comparative genomics databases and their data release. - More info: /documentation/info/comparas

data(*args, **kwargs)

GET info/data

Shows the data releases available on this REST server. May return more than one release (unfrequent non-standard Ensembl configuration). - More info: /documentation/info/data

eg_version(*args, **kwargs)

GET info/eg_version

Returns the Ensembl Genomes version of the databases backing this service - More info: /documentation/info/eg_version

external_dbs(*args, **kwargs)

GET info/external_dbs/:species

Lists all available external sources for a species. - More info: /documentation/info/external_dbs

family(*args, **kwargs)

GET family/id/:id

Retrieves gene family information by ID - More info: /documentation/info/family

family_member_id(*args, **kwargs)

GET family/member/id/:id

Retrieves gene families to which a gene belongs - More info: /documentation/info/family_member_id

family_member_symbol(*args, **kwargs)

GET family/member/symbol/:species/:symbol

Retrieves gene families to which a gene belongs - More info: /documentation/info/family_member_symbol

genetree(*args, **kwargs)

GET genetree/id/:id

Retrieves a gene tree dump for a gene tree stable identifier - More info: /documentation/info/genetree

genetree_member_id(*args, **kwargs)

GET genetree/member/id/:id

Retrieves a gene tree that contains the stable identifier - More info: /documentation/info/genetree_member_id

genetree_member_symbol(*args, **kwargs)

GET genetree/member/symbol/:species/:symbol

Retrieves a gene tree containing the gene identified by a symbol - More info: /documentation/info/genetree_member_symbol

genomic_alignment_region(*args, **kwargs)

GET alignment/region/:species/:region

Retrieves genomic alignments as separate blocks based on a region and species - More info: /documentation/info/genomic_alignment_region

homology_ensemblgene(*args, **kwargs)

GET homology/id/:id

Retrieves homology information (orthologs) by Ensembl gene id - More info: /documentation/info/homology_ensemblgene

homology_symbol(*args, **kwargs)

GET homology/symbol/:species/:symbol

Retrieves homology information (orthologs) by symbol - More info: /documentation/info/homology_symbol

info_divisions(*args, **kwargs)

GET info/divisions

Get list of all Ensembl divisions for which information is available - More info: /documentation/info/info_divisions

info_genome(*args, **kwargs)

GET info/genomes/:genome_name

Find information about a given genome - More info: /documentation/info/info_genome

info_genomes(*args, **kwargs)

GET info/genomes

Find information about all genomes. Response may be very large. - More info: /documentation/info/info_genomes

info_genomes_accession(*args, **kwargs)

GET info/genomes/accession/:accession

Find information about genomes containing a specified INSDC accession - More info: /documentation/info/info_genomes_accession

info_genomes_assembly(*args, **kwargs)

GET info/genomes/assembly/:assembly_id

Find information about a genome with a specified assembly - More info: /documentation/info/info_genomes_assembly

info_genomes_division(*args, **kwargs)

GET info/genomes/division/:division_name

Find information about all genomes in a given division. May be large for Ensembl Bacteria. - More info: /documentation/info/info_genomes_division

info_genomes_taxonomy(*args, **kwargs)

GET info/genomes/taxonomy/:taxon_name

Find information about all genomes beneath a given node of the taxonomy - More info: /documentation/info/info_genomes_taxonomy

list_all_microarrays(*args, **kwargs)

GET regulatory/species/:species/microarray

Returns information about all microarrays available for the given species - More info: /documentation/info/list_all_microarrays

lookup(*args, **kwargs)

GET lookup/id/:id

Find the species and database for a single identifier e.g. gene, transcript, protein - More info: /documentation/info/lookup

lookup_genome(*args, **kwargs)

GET lookup/genome/:name

Query for a named genome and retrieve the gene models - More info: /documentation/info/lookup_genome

lookup_post(*args, **kwargs)

POST lookup/id

Find the species and database for several identifiers. IDs that are not found are returned with no data. - More info: /documentation/info/lookup_post

ontology_ancestors(*args, **kwargs)

GET ontology/ancestors/:id

Reconstruct the entire ancestry of a term from is_a and part_of relationships - More info: /documentation/info/ontology_ancestors

ontology_ancestors_chart(*args, **kwargs)

GET ontology/ancestors/chart/:id

Reconstruct the entire ancestry of a term from is_a and part_of relationships. - More info: /documentation/info/ontology_ancestors_chart

ontology_descendants(*args, **kwargs)

GET ontology/descendants/:id

Find all the terms descended from a given term. By default searches are conducted within the namespace of the given identifier - More info: /documentation/info/ontology_descendants

ontology_id(*args, **kwargs)

GET ontology/id/:id

Search for an ontological term by its namespaced identifier - More info: /documentation/info/ontology_id

ontology_name(*args, **kwargs)

GET ontology/name/:name

Search for a list of ontological terms by their name - More info: /documentation/info/ontology_name

overlap_id(*args, **kwargs)

GET overlap/id/:id

Retrieves features (e.g. genes, transcripts, variants and more) that overlap a region defined by the given identifier. - More info: /documentation/info/overlap_id

overlap_region(*args, **kwargs)

GET overlap/region/:species/:region

Retrieves features (e.g. genes, transcripts, variants and more) that overlap a given region. - More info: /documentation/info/overlap_region

overlap_translation(*args, **kwargs)

GET overlap/translation/:id

Retrieve features related to a specific Translation as described by its stable ID (e.g. domains, variants). - More info: /documentation/info/overlap_translation

phenotype_gene(*args, **kwargs)

GET /phenotype/gene/:species/:gene

Return phenotype annotations for a given gene. - More info: /documentation/info/phenotype_gene

phenotype_region(*args, **kwargs)

GET /phenotype/region/:species/:region

Return phenotype annotations that overlap a given genomic region. - More info: /documentation/info/phenotype_region

ping(*args, **kwargs)

GET info/ping

Checks if the service is alive. - More info: /documentation/info/ping

probe_info(*args, **kwargs)

GET regulatory/species/:species/microarray/:microarray/probe/:probe

Returns information about a specific probe of a microarray - More info: /documentation/info/probe_info

rest(*args, **kwargs)

GET info/rest

Shows the current version of the Ensembl REST API. - More info: /documentation/info/rest

sequence_id(*args, **kwargs)

GET sequence/id/:id

Request multiple types of sequence by stable identifier. Supports feature masking and expand options. - More info: /documentation/info/sequence_id

sequence_id_post(*args, **kwargs)

POST sequence/id

Request multiple types of sequence by a stable identifier list. - More info: /documentation/info/sequence_id_post

sequence_region(*args, **kwargs)

GET sequence/region/:species/:region

Returns the genomic sequence of the specified region of the given species. Supports feature masking and expand options. - More info: /documentation/info/sequence_region

sequence_region_post(*args, **kwargs)

POST sequence/region/:species

Request multiple types of sequence by a list of regions. - More info: /documentation/info/sequence_region_post

software(*args, **kwargs)

GET info/software

Shows the current version of the Ensembl API used by the REST server. - More info: /documentation/info/software

species(*args, **kwargs)

GET info/species

Lists all available species, their aliases, available adaptor groups and data release. - More info: /documentation/info/species

symbol_lookup(*args, **kwargs)

GET lookup/symbol/:species/:symbol

Find the species and database for a symbol in a linked external database - More info: /documentation/info/symbol_lookup

symbol_post(*args, **kwargs)

POST lookup/symbol/:species/:symbol

Find the species and database for a set of symbols in a linked external database. Unknown symbols are omitted from the response. - More info: /documentation/info/symbol_post

taxonomy_classification(*args, **kwargs)

GET taxonomy/classification/:id

Return the taxonomic classification of a taxon node - More info: /documentation/info/taxonomy_classification

taxonomy_id(*args, **kwargs)

GET taxonomy/id/:id

Search for a taxonomic term by its identifier or name - More info: /documentation/info/taxonomy_id

taxonomy_name(*args, **kwargs)

GET taxonomy/name/:name

Search for a taxonomic id by a non-scientific name - More info: /documentation/info/taxonomy_name

variant_recoder(*args, **kwargs)

GET variant_recoder/:species/:id

Translate a variant identifier or HGVS notation to all possible variant IDs and HGVS - More info: /documentation/info/variant_recoder

variant_recoder_post(*args, **kwargs)

POST variant_recoder/:species

Translate a list of variant identifiers or HGVS notations to all possible variant IDs and HGVS - More info: /documentation/info/variant_recoder_post

variation(*args, **kwargs)

GET info/variation/:species

List the variation sources used in Ensembl for a species. - More info: /documentation/info/variation

variation_consequence_types(*args, **kwargs)

GET info/variation/consequence_types

Lists all variant consequence types. - More info: /documentation/info/variation_consequence_types

variation_id(*args, **kwargs)

GET variation/:species/:id

Uses a variant identifier (e.g. rsID) to return the variation features including optional genotype, phenotype and population data - More info: /documentation/info/variation_id

variation_populations(*args, **kwargs)

GET info/variation/populations/:species

List all populations for a species - More info: /documentation/info/variation_populations

variation_post(*args, **kwargs)

POST variation/:species/

Uses a list of variant identifiers (e.g. rsID) to return the variation features including optional genotype, phenotype and population data - More info: /documentation/info/variation_post

vep_id_get(*args, **kwargs)

GET vep/:species/id/:id

Fetch variant consequences based on a variant identifier - More info: /documentation/info/vep_id_get

vep_id_post(*args, **kwargs)

POST vep/:species/id

Fetch variant consequences for multiple ids - More info: /documentation/info/vep_id_post

vep_region_get(*args, **kwargs)

GET vep/:species/region/:region/:allele/

Fetch variant consequences - More info: /documentation/info/vep_region_get

vep_region_post(*args, **kwargs)

POST vep/:species/region

Fetch variant consequences for multiple regions - More info: /documentation/info/vep_region_post

xref_external(*args, **kwargs)

GET xrefs/symbol/:species/:symbol

Looks up an external symbol and returns all Ensembl objects linked to it. This can be a display name for a gene/transcript/translation, a synonym or an externally linked reference. If a gene’s transcript is linked to the supplied symbol the service will return both gene and transcript (it supports transient links). - More info: /documentation/info/xref_external

xref_id(*args, **kwargs)

GET xrefs/id/:id

Perform lookups of Ensembl Identifiers and retrieve their external references in other databases - More info: /documentation/info/xref_id

xref_name(*args, **kwargs)

GET xrefs/name/:species/:name

Performs a lookup based upon the primary accession or display label of an external reference and returning the information we hold about the entry - More info: /documentation/info/xref_name