InterMine is a powerful open source data warehouse system integrating diverse biological data sets (e.g. genomic, expression and protein data) for various organisms. Integrating data makes it possible to run sophisticated data mining queries that span domains of biological knowledge. A selected list of databases powered by InterMine is shown in Table 1:
Database | Organism | Data |
---|---|---|
FlyMine | Drosophila | Genes, homology, proteins, interactions, gene ontology, expression, regulation, phenotypes, pathways, diseases, resources, publications |
HumanMine | H. sapiens | Genomics, SNPs, GWAS, proteins, gene ontology, pathways, gene expression, interactions, publications, disease, orthologues, alleles |
MouseMine | M. musculus | Genomics, proteins, gene ontology, expression, interactions, pathways, phenotypes, diseases, homology, publications |
RatMine | R. norvegicus | Disease, gene ontology, genomics, interactions, phenotype, pathway, proteins, publication QTL, SNP |
WormMine | C. elegans | Genes, alleles, homology, go annotation, phenotypes, strains |
YeastMine | S. cerevisiae | Genomics, proteins, gene ontology, comparative genomics, phenotypes, interactions, literature, pathways, gene expression |
ZebrafishMine | D. rerio | Genes, constructs, disease, gene ontology, genotypes, homology, morpholinos, phenotypes |
TargetMine | H. sapiens, M. musculus | Genes, protein structures, chemical compounds, protein domains, gene function, pathways, interactions, disease, drug targets |
MitoMiner | H. sapiens, M. musculus, R. norvegicus, D. rerio, S. cerevisiae, S. pombe | Genes, homology, localisation evidence, Mitochondrial reference gene lists, phenotypes, diseases, expression, interactions, pathways, exome |
IndigoMine | Archae | Genomics |
ThaleMine | A. thaliana | Genomics, proteins, domains, homology, gene ontology, interactions, expression, publications, pathways, GeneRIF, stocks, phenotypes, alleles, insertions, TAIR |
MedicMine | Medicago truncatula | Genomics, pathways, gene ontology, publications, proteins, homology |
PhytoMine | over 50 plant genomes | Genes, proteins, expression, transcripts, homology |
Please see the InterMine home page for a full list of available InterMines.
InterMine includes an attractive, user-friendly web interface that works ‘out of the box’ and a powerful, scriptable web-service API to allow programmatic access to your data. This R package provides an interface with the InterMine-powered databases through Web services.
Let’s start with a simple task - find the homologues of gene ABO.
First, we look at what databases are available.
library(InterMineR)
listMines()
## BMAP
## "https://bmap.jgi.doe.gov/bmapmine/"
## BeanMine
## "https://mines.legumeinfo.org/beanmine"
## BovineMine
## "http://bovinegenome.org/bovinemine"
## CHOmine
## "https://chomine.boku.ac.at/chomine"
## ChickpeaMine
## "https://mines.legumeinfo.org/chickpeamine"
## CowpeaMine
## "https://mines.legumeinfo.org/cowpeamine"
## FlyMine
## "http://www.flymine.org/flymine"
## GrapeMine
## "http://urgi.versailles.inra.fr/GrapeMine"
## HumanMine
## "http://www.humanmine.org/humanmine"
## HymenopteraMine
## "http://hymenopteragenome.org/hymenopteramine"
## IndigoMine
## "http://www.cbrc.kaust.edu.sa/indigo"
## LegumeMine
## "https://mines.legumeinfo.org/legumemine"
## MaizeMine
## "http://maizemine.rnet.missouri.edu:8080/maizemine"
## MedicMine
## "http://medicmine.jcvi.org/medicmine"
## MitoMiner
## "http://mitominer.mrc-mbu.cam.ac.uk/release-4.0"
## ModMine
## "http://intermine.modencode.org/release-33"
## MouseMine
## "http://www.mousemine.org/mousemine"
## OakMine
## "https://urgi.versailles.inra.fr/OakMine_PM1N"
## PeanutMine
## "https://mines.legumeinfo.org/peanutmine"
## PhytoMine
## "https://phytozome.jgi.doe.gov/phytomine/"
## PlanMine
## "http://planmine.mpi-cbg.de/planmine"
## RatMine
## "http://ratmine.mcw.edu/ratmine"
## RepetDB
## "http://urgi.versailles.inra.fr/repetdb"
## SoyMine
## "https://mines.legumeinfo.org/soymine"
## TargetMine
## "http://targetmine.mizuguchilab.org/targetmine"
## TetraMine
## "http://adenine.bradley.edu/tetramine"
## ThaleMine
## "https://apps.araport.org/thalemine"
## WheatMine
## "https://urgi.versailles.inra.fr/WheatMine"
## WormMine
## "http://intermine.wormbase.org/tools/wormmine/"
## XenMine
## "http://www.xenmine.org/xenmine"
## YeastMine
## "https://yeastmine.yeastgenome.org/yeastmine"
## ZebrafishMine
## "http://zebrafishmine.org"
Since we would like to query human genes, we select HumanMine.
# load HumaMine
im <- initInterMine(mine=listMines()["HumanMine"])
im
## $mine
## HumanMine
## "http://www.humanmine.org/humanmine"
##
## $token
## [1] ""
Both in InterMine database website and in InterMineR, you are able to build custom queries. However, to facilitate the retrieval of information from InterMine databases, a variety of pre-built queries, called templates, have also been made available. Templates are queries that have already been created with a fixed set of output columns and one or more constraints.
# Get template (collection of pre-defined queries)
template = getTemplates(im)
head(template)
## name
## 1 Tissue_Expression_illumina
## 2 humDisGeneOrthol2
## 3 PhenotypeGene
## 4 Protein_complex_details
## 5 disExprGene
## 6 protein_interactions
## title
## 1 Tissue --> Gene Expression (Illumina body map)
## 2 Human Disease --> Human Gene + Orthologue Gene(s)
## 3 Mouse Phenotype --> Mouse Genes + Orthologous genes
## 4 Protein --> Protein Complex
## 5 Disease Expression --> Genes
## 6 Protein --> Interactions
We would like to find templates involving genes.
# Get gene-related templates
template[grep("gene", template$name, ignore.case=TRUE),]
## name
## 2 humDisGeneOrthol2
## 3 PhenotypeGene
## 5 disExprGene
## 7 Gene_Interactions2
## 8 Protein_Gene_Ortho
## 9 GOterm_Gene
## 11 Gene_Alleles_Disease2
## 13 ChromRegion_GenesTransExon
## 16 GeneExpress
## 17 Disease_Genes2
## 18 Gene_Location
## 19 Protein_GeneChromosomeLength
## 20 Gene_Identifiers
## 22 Gene_Pathway
## 23 gene_complex_details
## 26 PathwayGenes
## 28 Gene_Protein
## 29 Gene_OverlapppingGenes
## 32 Gene_To_Publications
## 33 Gene_Interactions_forReportPage
## 35 Gene_Disease_HPO
## 36 Gene_GO
## 38 GeneInteractorsExpression
## 39 Gene_particularGoannotation
## 40 Gene_TissueExpressionIllumina
## 41 Gene_HPOphenotype_2
## 42 domain_protein_gene
## 43 Gene_Expression_GTex
## 44 Pathway_ProteinGene
## 47 GeneHPOparent_Genes_2
## 48 Gene_proteindomain
## 49 Gene_inGWAS
## 50 geneGWAS_reportPg
## 51 geneInteractiongene
## 52 Gene_Disease2
## 53 Term_inGWASoptionalGene
## 54 Gene_proteinAtlasExpression2
## 55 GeneOrthAllele
## 59 Gene_Orth
## 60 ChromRegion_Genes
## 62 GenePathway_interactions2
## 63 Gene_AllelePhen
## title
## 2 Human Disease --> Human Gene + Orthologue Gene(s)
## 3 Mouse Phenotype --> Mouse Genes + Orthologous genes
## 5 Disease Expression --> Genes
## 7 Gene --> Interactions
## 8 Protein --> Gene and Orthologues
## 9 GO term --> Genes
## 11 Gene --> Alleles and Disease
## 13 Chromosomal Location --> All Genes + Transcripts + Exons
## 16 Gene --> Gene Expression (Tissue, Disease; Array Express, E-MTAB-62)
## 17 Disease --> Gene(s)
## 18 Gene --> Chromosomal location.
## 19 Protein --> Gene.
## 20 Gene --> All identifiers.
## 22 Gene --> Pathway
## 23 Gene --> protein complex
## 26 Pathway --> Genes
## 28 Gene --> Proteins.
## 29 Gene --> Overlapping genes.
## 32 Gene --> Publications.
## 33 Gene --> Physical and Genetic Interactions
## 35 Gene --> Disease + HPO annotations (Human Phenotype Ontology
## 36 Gene --> GO terms.
## 38 Gene + Tissue Expression --> Interactors that are expressed in that tissue
## 39 Gene + GO term --> Genes by GO term
## 40 Gene --> Tissue Expression (Illumina body map)
## 41 Gene -> HPO annotation (Human Phenotype Ontology)
## 42 Protein Domain --> Protein and Genes
## 43 Gene --> Tissue Expression (GTex data)
## 44 Pathway --> Protein and Gene
## 47 Gene + HPO Phenotype parent term -> Genes
## 48 Gene --> Protein + Domains
## 49 Gene --> GWAS hit
## 50 Gene Report --> GWAS hit
## 51 Gene A --> Interaction <-- Gene B
## 52 Gene --> Disease (OMIM)
## 53 GWAS term --> SNP + Associated gene if available
## 54 Gene --> Protein tissue Localisation
## 55 Gene (Hum OR Rat) --> Mouse Allele (Phenotype)
## 59 Gene --> Orthologues
## 60 Region --> Genes
## 62 Gene + Pathway --> Interactions
## 63 Mouse Gene --> Allele [Phenotype]
The template Gene_Orth seems to be what we want. Let’s look at this template in more detail.
# Query for gene orthologs
queryGeneOrth = getTemplateQuery(
im = im,
name = "Gene_Orth"
)
queryGeneOrth
## $model
## name
## "genomic"
##
## $title
## [1] "Gene --> Orthologues"
##
## $description
## [1] "For a given Gene (or List of Genes) in named organism (default: Human) returns the orthologues in a different organisms. [keywords: homologue, homolog, paralogue, paralogue, ortholog]"
##
## $select
## [1] "Gene.primaryIdentifier"
## [2] "Gene.symbol"
## [3] "Gene.homologues.homologue.primaryIdentifier"
## [4] "Gene.homologues.homologue.symbol"
## [5] "Gene.homologues.homologue.organism.shortName"
##
## $name
## [1] "Gene_Orth"
##
## $comment
## [1] ""
##
## $orderBy
## $orderBy[[1]]
## Gene.symbol
## "ASC"
##
##
## $where
## $where[[1]]
## $where[[1]]$path
## [1] "Gene"
##
## $where[[1]]$op
## [1] "LOOKUP"
##
## $where[[1]]$code
## [1] "A"
##
## $where[[1]]$editable
## [1] TRUE
##
## $where[[1]]$switchable
## [1] FALSE
##
## $where[[1]]$switched
## [1] "LOCKED"
##
## $where[[1]]$value
## [1] "PPARG"
##
## $where[[1]]$extraValue
## [1] "H. sapiens"
There are three essential members in a query - SELECT, WHERE and constraintLogic.
What does ‘Gene.symbol’ mean? What is ‘Gene.homologues.homologue.symbol’?
Let’s take a look at the data model.
model <- getModel(im)
head(model)
## type child_name child_type
## 1 Allele Alternate
## 2 Allele Clinical Significance
## 3 Allele Id
## 4 Allele Name
## 5 Allele Primary Identifier
## 6 Allele Reference
Let’s look at the children of the Gene data type.
model[which(model$type=="Gene"),]
## type child_name child_type
## 924 Gene Brief Description
## 925 Gene Cytological Location
## 926 Gene Description
## 927 Gene Id
## 928 Gene Length
## 929 Gene Name
## 930 Gene Primary Identifier
## 931 Gene Score
## 932 Gene Score Type
## 933 Gene Secondary Identifier
## 934 Gene Symbol
## 935 Gene alleles Allele
## 936 Gene atlasExpression AtlasExpression
## 937 Gene CDSs CDS
## 938 Gene chromosome Chromosome
## 939 Gene crossReferences CrossReference
## 940 Gene dataSets DataSet
## 941 Gene diseases Disease
## 942 Gene exons Exon
## 943 Gene goAnnotation GOAnnotation
## 944 Gene flankingRegions GeneFlankingRegion
## 945 Gene homologues Homologue
## 946 Gene interactions Interaction
## 947 Gene downstreamIntergenicRegion IntergenicRegion
## 948 Gene upstreamIntergenicRegion IntergenicRegion
## 949 Gene introns Intron
## 950 Gene chromosomeLocation Location
## 951 Gene locatedFeatures Location
## 952 Gene locations Location
## 953 Gene ontologyAnnotations OntologyAnnotation
## 954 Gene organism Organism
## 955 Gene pathways Pathway
## 956 Gene probeSets ProbeSet
## 957 Gene proteins Protein
## 958 Gene proteinAtlasExpression ProteinAtlasExpression
## 959 Gene publications Publication
## 960 Gene rnaSeqResults RNASeqResult
## 961 Gene regulatoryRegions RegulatoryRegion
## 962 Gene SNPs SNP
## 963 Gene sequenceOntologyTerm SOTerm
## 964 Gene sequence Sequence
## 965 Gene childFeatures SequenceFeature
## 966 Gene overlappingFeatures SequenceFeature
## 967 Gene strain Strain
## 968 Gene synonyms Synonym
## 969 Gene transcripts Transcript
## 970 Gene UTRs UTR
Gene has a field called ‘symbol’ (hence the column Gene.symbol). Gene also has a child called homologues, which is of the Homologue data type.
model[which(model$type=="Homologue"),]
## type child_name child_type
## 1091 Homologue Id
## 1092 Homologue Type
## 1093 Homologue crossReferences CrossReference
## 1094 Homologue dataSets DataSet
## 1095 Homologue gene Gene
## 1096 Homologue homologue Gene
## 1097 Homologue evidence OrthologueEvidence
Homologue has a child called ‘gene’ which is of the type ‘Gene’, which we saw above has a field called ‘symbol’ (hence the column Gene.homologues.homologue.symbol).
Let’s now run our template.
resGeneOrth <- runQuery(im, queryGeneOrth)
head(resGeneOrth)
## Gene.primaryIdentifier Gene.symbol
## 1 5468 PPARG
## 2 5468 PPARG
## 3 5468 PPARG
## 4 5468 PPARG
## 5 5468 PPARG
## 6 5468 PPARG
## Gene.homologues.homologue.primaryIdentifier
## 1 10062
## 2 5465
## 3 5467
## 4 5914
## 5 5915
## 6 5916
## Gene.homologues.homologue.symbol
## 1 NR1H3
## 2 PPARA
## 3 PPARD
## 4 RARA
## 5 RARB
## 6 RARG
## Gene.homologues.homologue.organism.shortName
## 1 H. sapiens
## 2 H. sapiens
## 3 H. sapiens
## 4 H. sapiens
## 5 H. sapiens
## 6 H. sapiens
Let’s modify the query to find the orthologues of the gene ABO. We want to change the ‘value’ attribute from PPARG to ABO.
There are two ways to build a query in InterMineR.
We can either build a query as a list object with newQuery
function, and assign all input values (selection of retrieved data type, constraints, etc.) as items of that list,
Or we can build the query as an InterMineR-class
object with the functions setConstraint
, which allows us to generate a new or modify an existing list of constraints, and setQuery
, which generates the query as a InterMineR-class
object.
setConstraints
and setQuery
functions are designed to facilitate the generation of queries for InterMine instances and avoid using multiple iterative loops, especially when it is required to include multiple constraints or constraint values (e.g. genes, organisms) in your query.
# modify directly the value of the first constraint from the list query
queryGeneOrth$where[[1]][["value"]] <- "ABO"
# or modify the value of the first constraint from the list query with setConstraints
queryGeneOrth$where = setConstraints(
modifyQueryConstraints = queryGeneOrth,
m.index = 1,
values = list("ABO")
)
queryGeneOrth$where
## [[1]]
## [[1]]$path
## [1] "Gene"
##
## [[1]]$op
## [1] "LOOKUP"
##
## [[1]]$code
## [1] "A"
##
## [[1]]$editable
## [1] TRUE
##
## [[1]]$switchable
## [1] FALSE
##
## [[1]]$switched
## [1] "LOCKED"
##
## [[1]]$value
## [1] "ABO"
##
## [[1]]$extraValue
## [1] "H. sapiens"
Note the value is now equal to ‘ABO’. Let’s rerun our query with the new constraint.
resGeneOrth <- runQuery(im, queryGeneOrth)
head(resGeneOrth)
## Gene.primaryIdentifier Gene.symbol
## 1 28 ABO
## 2 28 ABO
## 3 28 ABO
## 4 28 ABO
## 5 28 ABO
## 6 28 ABO
## Gene.homologues.homologue.primaryIdentifier
## 1 127550
## 2 26301
## 3 360203
## 4 MGI:2135738
## 5 RGD:1311152
## 6 RGD:2307241
## Gene.homologues.homologue.symbol
## 1 A3GALT2
## 2 GBGT1
## 3 GLT6D1
## 4 Abo
## 5 Abo2
## 6 Abo
## Gene.homologues.homologue.organism.shortName
## 1 H. sapiens
## 2 H. sapiens
## 3 H. sapiens
## 4 M. musculus
## 5 R. norvegicus
## 6 R. norvegicus
Now we are seeing orthologues for the ABO gene. Let’s add the organism to the view to make sure we are looking at the desired gene.
You can also add additional filters. Let’s exclude all homologues where organism is H. sapiens.
There are four parts of a constraint to add:
newConstraint <- list(
path=c("Gene.homologues.homologue.organism.shortName"),
op=c("!="),
value=c("H. sapiens"),
code=c("B")
)
queryGeneOrth$where[[2]] <- newConstraint
queryGeneOrth$where
## [[1]]
## [[1]]$path
## [1] "Gene"
##
## [[1]]$op
## [1] "LOOKUP"
##
## [[1]]$code
## [1] "A"
##
## [[1]]$editable
## [1] TRUE
##
## [[1]]$switchable
## [1] FALSE
##
## [[1]]$switched
## [1] "LOCKED"
##
## [[1]]$value
## [1] "ABO"
##
## [[1]]$extraValue
## [1] "H. sapiens"
##
##
## [[2]]
## [[2]]$path
## [1] "Gene.homologues.homologue.organism.shortName"
##
## [[2]]$op
## [1] "!="
##
## [[2]]$value
## [1] "H. sapiens"
##
## [[2]]$code
## [1] "B"
Our new filter has been added successfully. Rerun the query and you see you only have non-Homo sapiens orthologues.
resGeneOrth <- runQuery(im, queryGeneOrth)
resGeneOrth
## Gene.primaryIdentifier Gene.symbol
## 1 28 ABO
## 2 28 ABO
## 3 28 ABO
## 4 28 ABO
## 5 28 ABO
## 6 28 ABO
## 7 28 ABO
## 8 28 ABO
## 9 28 ABO
## 10 28 ABO
## 11 28 ABO
## 12 28 ABO
## Gene.homologues.homologue.primaryIdentifier
## 1 MGI:2135738
## 2 RGD:1311152
## 3 RGD:2307241
## 4 RGD:628609
## 5 ZDB-GENE-031204-4
## 6 ZDB-GENE-040426-1117
## 7 ZDB-GENE-040912-46
## 8 ZDB-GENE-060531-15
## 9 ZDB-GENE-060531-59
## 10 ZDB-GENE-060531-70
## 11 ZDB-GENE-060531-71
## 12 ZDB-GENE-081104-23
## Gene.homologues.homologue.symbol
## 1 Abo
## 2 Abo2
## 3 Abo
## 4 Abo3
## 5
## 6
## 7
## 8
## 9
## 10
## 11
## 12
## Gene.homologues.homologue.organism.shortName
## 1 M. musculus
## 2 R. norvegicus
## 3 R. norvegicus
## 4 R. norvegicus
## 5 D. rerio
## 6 D. rerio
## 7 D. rerio
## 8 D. rerio
## 9 D. rerio
## 10 D. rerio
## 11 D. rerio
## 12 D. rerio
You can also add additional columns to the output. For instance, where do these homologues come from? Let’s add this information.
Let’s see what we know about homologues.
model[which(model$type=="Homologue"),]
## type child_name child_type
## 1091 Homologue Id
## 1092 Homologue Type
## 1093 Homologue crossReferences CrossReference
## 1094 Homologue dataSets DataSet
## 1095 Homologue gene Gene
## 1096 Homologue homologue Gene
## 1097 Homologue evidence OrthologueEvidence
The Homologue data type has an ‘dataSets’ reference of type ‘DataSet’.
model[which(model$type=="DataSet"),]
## type child_name child_type
## 648 DataSet Description
## 649 DataSet Id
## 650 DataSet Name
## 651 DataSet URL
## 652 DataSet Version
## 653 DataSet bioEntities BioEntity
## 654 DataSet dataSource DataSource
## 655 DataSet publication Publication
DataSet has a child called name. Add Gene.homologues.dataSets.name to the view. We’ll add it as the last column, we can see from above there are 5 other columns already so we’ll put it as #6:
# use setQuery function which will create an InterMineR-class query
queryGeneOrth.InterMineR = setQuery(
inheritQuery = queryGeneOrth,
select = c(queryGeneOrth$select,
"Gene.homologues.dataSets.name")
)
getSelect(queryGeneOrth.InterMineR)
## [1] "Gene.primaryIdentifier"
## [2] "Gene.symbol"
## [3] "Gene.homologues.homologue.primaryIdentifier"
## [4] "Gene.homologues.homologue.symbol"
## [5] "Gene.homologues.homologue.organism.shortName"
## [6] "Gene.homologues.dataSets.name"
#queryGeneOrth.InterMineR@select
# or assign new column directly to the existing list query
queryGeneOrth$select[[6]] <- "Gene.homologues.dataSets.name"
queryGeneOrth$select
## [1] "Gene.primaryIdentifier"
## [2] "Gene.symbol"
## [3] "Gene.homologues.homologue.primaryIdentifier"
## [4] "Gene.homologues.homologue.symbol"
## [5] "Gene.homologues.homologue.organism.shortName"
## [6] "Gene.homologues.dataSets.name"
# run queries
resGeneOrth.InterMineR <- runQuery(im, queryGeneOrth.InterMineR)
resGeneOrth <- runQuery(im, queryGeneOrth)
all(resGeneOrth == resGeneOrth.InterMineR)
## [1] TRUE
head(resGeneOrth, 3)
## Gene.primaryIdentifier Gene.symbol
## 1 28 ABO
## 2 28 ABO
## 3 28 ABO
## Gene.homologues.homologue.primaryIdentifier
## 1 MGI:2135738
## 2 RGD:1311152
## 3 RGD:2307241
## Gene.homologues.homologue.symbol
## 1 Abo
## 2 Abo2
## 3 Abo
## Gene.homologues.homologue.organism.shortName
## 1 M. musculus
## 2 R. norvegicus
## 3 R. norvegicus
## Gene.homologues.dataSets.name
## 1 Orthologue and paralogue predictions
## 2 Orthologue and paralogue predictions
## 3 Orthologue and paralogue predictions
NB: adding columns can result in changing the row count.
The constraintLogic, if not given, is ‘A and B’. We would now try to explicitly specify the constraintLogic. A and B corresponds to the ‘code’ for each constraint.
queryGeneOrth$constraintLogic <- "A and B"
queryGeneOrth$constraintLogic
## [1] "A and B"
Run the query again and see no change:
resGeneOrth <- runQuery(im, queryGeneOrth)
resGeneOrth
## Gene.primaryIdentifier Gene.symbol
## 1 28 ABO
## 2 28 ABO
## 3 28 ABO
## 4 28 ABO
## 5 28 ABO
## 6 28 ABO
## 7 28 ABO
## 8 28 ABO
## 9 28 ABO
## 10 28 ABO
## 11 28 ABO
## 12 28 ABO
## Gene.homologues.homologue.primaryIdentifier
## 1 MGI:2135738
## 2 RGD:1311152
## 3 RGD:2307241
## 4 RGD:628609
## 5 ZDB-GENE-031204-4
## 6 ZDB-GENE-040426-1117
## 7 ZDB-GENE-040912-46
## 8 ZDB-GENE-060531-15
## 9 ZDB-GENE-060531-59
## 10 ZDB-GENE-060531-70
## 11 ZDB-GENE-060531-71
## 12 ZDB-GENE-081104-23
## Gene.homologues.homologue.symbol
## 1 Abo
## 2 Abo2
## 3 Abo
## 4 Abo3
## 5
## 6
## 7
## 8
## 9
## 10
## 11
## 12
## Gene.homologues.homologue.organism.shortName
## 1 M. musculus
## 2 R. norvegicus
## 3 R. norvegicus
## 4 R. norvegicus
## 5 D. rerio
## 6 D. rerio
## 7 D. rerio
## 8 D. rerio
## 9 D. rerio
## 10 D. rerio
## 11 D. rerio
## 12 D. rerio
## Gene.homologues.dataSets.name
## 1 Orthologue and paralogue predictions
## 2 Orthologue and paralogue predictions
## 3 Orthologue and paralogue predictions
## 4 Orthologue and paralogue predictions
## 5 Orthologue and paralogue predictions
## 6 Orthologue and paralogue predictions
## 7 Orthologue and paralogue predictions
## 8 Orthologue and paralogue predictions
## 9 Orthologue and paralogue predictions
## 10 Orthologue and paralogue predictions
## 11 Orthologue and paralogue predictions
## 12 Orthologue and paralogue predictions
Change to be ‘A or B’ and see how the results change.
- Start with the template Gene GO
queryGeneGO <- getTemplateQuery(im, "Gene_GO")
queryGeneGO
## $model
## name
## "genomic"
##
## $title
## [1] "Gene --> GO terms."
##
## $description
## [1] "Search for GO annotations for a particular gene (or List of Genes)."
##
## $select
## [1] "Gene.primaryIdentifier"
## [2] "Gene.symbol"
## [3] "Gene.goAnnotation.ontologyTerm.identifier"
## [4] "Gene.goAnnotation.ontologyTerm.name"
## [5] "Gene.goAnnotation.ontologyTerm.namespace"
## [6] "Gene.goAnnotation.evidence.code.code"
## [7] "Gene.goAnnotation.ontologyTerm.parents.identifier"
## [8] "Gene.goAnnotation.ontologyTerm.parents.name"
## [9] "Gene.goAnnotation.qualifier"
##
## $name
## [1] "Gene_GO"
##
## $comment
## [1] "Added 15NOV2010: ML"
##
## $orderBy
## $orderBy[[1]]
## Gene.primaryIdentifier
## "ASC"
##
##
## $where
## $where[[1]]
## $where[[1]]$path
## [1] "Gene"
##
## $where[[1]]$op
## [1] "LOOKUP"
##
## $where[[1]]$code
## [1] "A"
##
## $where[[1]]$editable
## [1] TRUE
##
## $where[[1]]$switchable
## [1] FALSE
##
## $where[[1]]$switched
## [1] "LOCKED"
##
## $where[[1]]$value
## [1] "PPARG"
##
## $where[[1]]$extraValue
## [1] "H. sapiens"
- Modify the view to display a compact view
queryGeneGO$select <- queryGeneGO$select[2:5]
queryGeneGO$select
## [1] "Gene.symbol"
## [2] "Gene.goAnnotation.ontologyTerm.identifier"
## [3] "Gene.goAnnotation.ontologyTerm.name"
## [4] "Gene.goAnnotation.ontologyTerm.namespace"
- Modify the constraints to look for gene ABO.
queryGeneGO$where[[1]][["value"]] <- "ABO"
queryGeneGO$where
## [[1]]
## [[1]]$path
## [1] "Gene"
##
## [[1]]$op
## [1] "LOOKUP"
##
## [[1]]$code
## [1] "A"
##
## [[1]]$editable
## [1] TRUE
##
## [[1]]$switchable
## [1] FALSE
##
## [[1]]$switched
## [1] "LOCKED"
##
## [[1]]$value
## [1] "ABO"
##
## [[1]]$extraValue
## [1] "H. sapiens"
- Run the query
resGeneGO <- runQuery(im, queryGeneGO )
head(resGeneGO)
## Gene.symbol Gene.goAnnotation.ontologyTerm.identifier
## 1 ABO GO:0000166
## 2 ABO GO:0003823
## 3 ABO GO:0004380
## 4 ABO GO:0004381
## 5 ABO GO:0005576
## 6 ABO GO:0005794
## Gene.goAnnotation.ontologyTerm.name
## 1 nucleotide binding
## 2 antigen binding
## 3 glycoprotein-fucosylgalactoside alpha-N-acetylgalactosaminyltransferase activity
## 4 fucosylgalactoside 3-alpha-galactosyltransferase activity
## 5 extracellular region
## 6 Golgi apparatus
## Gene.goAnnotation.ontologyTerm.namespace
## 1 molecular_function
## 2 molecular_function
## 3 molecular_function
## 4 molecular_function
## 5 cellular_component
## 6 cellular_component
- Start with the template Gene GO
queryGOGene <- getTemplateQuery(im, "GOterm_Gene")
queryGOGene
## $model
## name
## "genomic"
##
## $title
## [1] "GO term --> Genes"
##
## $description
## [1] "Search for Genes in a specified organism that are associated with a particular Gene Ontology (GO) annotation."
##
## $select
## [1] "Gene.primaryIdentifier"
## [2] "Gene.symbol"
## [3] "Gene.name"
## [4] "Gene.goAnnotation.ontologyTerm.identifier"
## [5] "Gene.goAnnotation.ontologyTerm.name"
## [6] "Gene.organism.shortName"
##
## $constraintLogic
## [1] "A and B"
##
## $name
## [1] "GOterm_Gene"
##
## $comment
## [1] "Added 26OCT2010: ML"
##
## $orderBy
## $orderBy[[1]]
## Gene.symbol
## "ASC"
##
##
## $where
## $where[[1]]
## $where[[1]]$path
## [1] "Gene.goAnnotation.ontologyTerm.name"
##
## $where[[1]]$op
## [1] "LIKE"
##
## $where[[1]]$code
## [1] "A"
##
## $where[[1]]$editable
## [1] TRUE
##
## $where[[1]]$switchable
## [1] FALSE
##
## $where[[1]]$switched
## [1] "LOCKED"
##
## $where[[1]]$value
## [1] "DNA binding"
##
##
## $where[[2]]
## $where[[2]]$path
## [1] "Gene.organism.shortName"
##
## $where[[2]]$op
## [1] "="
##
## $where[[2]]$code
## [1] "B"
##
## $where[[2]]$editable
## [1] TRUE
##
## $where[[2]]$switchable
## [1] FALSE
##
## $where[[2]]$switched
## [1] "LOCKED"
##
## $where[[2]]$value
## [1] "H. sapiens"
- Modify the view to display a compact view
queryGOGene$select <- queryGOGene$select[2:5]
queryGOGene$select
## [1] "Gene.symbol"
## [2] "Gene.name"
## [3] "Gene.goAnnotation.ontologyTerm.identifier"
## [4] "Gene.goAnnotation.ontologyTerm.name"
- Modify the constraints to look for GO term ‘metal ion binding’
queryGOGene$where[[1]]$value = "metal ion binding"
queryGOGene$where
## [[1]]
## [[1]]$path
## [1] "Gene.goAnnotation.ontologyTerm.name"
##
## [[1]]$op
## [1] "LIKE"
##
## [[1]]$code
## [1] "A"
##
## [[1]]$editable
## [1] TRUE
##
## [[1]]$switchable
## [1] FALSE
##
## [[1]]$switched
## [1] "LOCKED"
##
## [[1]]$value
## [1] "metal ion binding"
##
##
## [[2]]
## [[2]]$path
## [1] "Gene.organism.shortName"
##
## [[2]]$op
## [1] "="
##
## [[2]]$code
## [1] "B"
##
## [[2]]$editable
## [1] TRUE
##
## [[2]]$switchable
## [1] FALSE
##
## [[2]]$switched
## [1] "LOCKED"
##
## [[2]]$value
## [1] "H. sapiens"
- Run the query
resGOGene <- runQuery(im, queryGOGene )
head(resGOGene)
## Gene.symbol Gene.name
## 1 A3GALT2 alpha 1,3-galactosyltransferase 2
## 2 AARS alanyl-tRNA synthetase
## 3 AARS2 alanyl-tRNA synthetase 2, mitochondrial
## 4 AARSD1 alanyl-tRNA synthetase domain containing 1
## 5 ABAT 4-aminobutyrate aminotransferase
## 6 ABCG5 ATP binding cassette subfamily G member 5
## Gene.goAnnotation.ontologyTerm.identifier
## 1 GO:0046872
## 2 GO:0046872
## 3 GO:0046872
## 4 GO:0046872
## 5 GO:0046872
## 6 GO:0046872
## Gene.goAnnotation.ontologyTerm.name
## 1 metal ion binding
## 2 metal ion binding
## 3 metal ion binding
## 4 metal ion binding
## 5 metal ion binding
## 6 metal ion binding
- Start with the Gene_Location template, update to search for ABCA6 gene.
queryGeneLoc = getTemplateQuery(im, "Gene_Location")
queryGeneLoc$where[[2]][["value"]] = "ABCA6"
resGeneLoc= runQuery(im, queryGeneLoc)
resGeneLoc
## Gene.primaryIdentifier Gene.secondaryIdentifier Gene.symbol
## 1 23460 ENSG00000154262 ABCA6
## Gene.name
## 1 ATP binding cassette subfamily A member 6
## Gene.chromosome.primaryIdentifier Gene.locations.start Gene.locations.end
## 1 17 69062044 69141927
## Gene.locations.strand
## 1 -1
We’re going to use the output (gene location) as input for the next query.
- Define a new query
# set constraints
constraints = setConstraints(
paths = c(
"Gene.chromosome.primaryIdentifier",
"Gene.locations.start",
"Gene.locations.end",
"Gene.organism.name"
),
operators = c(
"=",
">=",
"<=",
"="
),
values = list(
resGeneLoc[1, "Gene.chromosome.primaryIdentifier"],
as.character(as.numeric(resGeneLoc[1, "Gene.locations.start"])-50000),
as.character(as.numeric(resGeneLoc[1, "Gene.locations.end"])+50000),
"Homo sapiens"
)
)
# set InterMineR-class query
queryNeighborGene = setQuery(
select = c("Gene.primaryIdentifier",
"Gene.symbol",
"Gene.chromosome.primaryIdentifier",
"Gene.locations.start",
"Gene.locations.end",
"Gene.locations.strand"),
where = constraints
)
summary(queryNeighborGene)
## path op value code
## 1 Gene.chromosome.primaryIdentifier = 17 A
## 2 Gene.locations.start >= 69012044 B
## 3 Gene.locations.end <= 69191927 C
## 4 Gene.organism.name = Homo sapiens D
- Run the query
resNeighborGene <- runQuery(im, queryNeighborGene)
resNeighborGene
## Gene.primaryIdentifier Gene.symbol Gene.chromosome.primaryIdentifier
## 1 100616316 MIR4524A 17
## 2 100847008 MIR4524B 17
## 3 23460 ABCA6 17
## Gene.locations.start Gene.locations.end Gene.locations.strand
## 1 69099564 69099632 -1
## 2 69099542 69099656 1
## 3 69062044 69141927 -1
- Plot the genes
resNeighborGene$Gene.locations.strand[which(resNeighborGene$Gene.locations.strand==1)]="+"
resNeighborGene$Gene.locations.strand[which(resNeighborGene$Gene.locations.strand==-1)]="-"
gene.idx = which(nchar(resNeighborGene$Gene.symbol)==0)
resNeighborGene$Gene.symbol[gene.idx]=resNeighborGene$Gene.primaryIdentifier[gene.idx]
require(Gviz)
annTrack = AnnotationTrack(
start=resNeighborGene$Gene.locations.start,
end=resNeighborGene$Gene.locations.end,
strand=resNeighborGene$Gene.locations.strand,
chromosome=resNeighborGene$Gene.chromosome.primaryIdentifier[1],
genome="GRCh38",
name="around ABCA6",
id=resNeighborGene$Gene.symbol)
gtr <- GenomeAxisTrack()
itr <- IdeogramTrack(genome="hg38", chromosome="chr17")
plotTracks(list(gtr, itr, annTrack), shape="box", showFeatureId=TRUE, fontcolor="black")
sessionInfo()
## R version 3.5.1 Patched (2018-07-12 r74967)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 16.04.5 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.8-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.8-bioc/R/lib/libRlapack.so
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=C LC_COLLATE=C
## [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
## [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
## [9] LC_ADDRESS=C LC_TELEPHONE=C
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
##
## attached base packages:
## [1] grid parallel stats4 stats graphics grDevices utils
## [8] datasets methods base
##
## other attached packages:
## [1] Gviz_1.26.3 GenomicRanges_1.34.0 GenomeInfoDb_1.18.1
## [4] org.Hs.eg.db_3.7.0 GO.db_3.7.0 GeneAnswers_2.24.0
## [7] RColorBrewer_1.1-2 Heatplus_2.28.0 MASS_7.3-51.1
## [10] annotate_1.60.0 XML_3.98-1.16 AnnotationDbi_1.44.0
## [13] IRanges_2.16.0 S4Vectors_0.20.1 Biobase_2.42.0
## [16] BiocGenerics_0.28.0 RCurl_1.95-4.11 bitops_1.0-6
## [19] igraph_1.2.2 RSQLite_2.1.1 InterMineR_1.4.1
## [22] BiocStyle_2.10.0
##
## loaded via a namespace (and not attached):
## [1] colorspace_1.3-2 rprojroot_1.3-2
## [3] biovizBase_1.30.0 htmlTable_1.12
## [5] XVector_0.22.0 base64enc_0.1-3
## [7] dichromat_2.0-0 rstudioapi_0.8
## [9] bit64_0.9-7 sqldf_0.4-11
## [11] xml2_1.2.0 splines_3.5.1
## [13] knitr_1.20 Formula_1.2-3
## [15] jsonlite_1.5 Rsamtools_1.34.0
## [17] cluster_2.0.7-1 graph_1.60.0
## [19] BiocManager_1.30.4 compiler_3.5.1
## [21] httr_1.3.1 backports_1.1.2
## [23] assertthat_0.2.0 Matrix_1.2-15
## [25] lazyeval_0.2.1 acepack_1.4.1
## [27] htmltools_0.3.6 prettyunits_1.0.2
## [29] tools_3.5.1 bindrcpp_0.2.2
## [31] gtable_0.2.0 glue_1.3.0
## [33] GenomeInfoDbData_1.2.0 dplyr_0.7.8
## [35] Rcpp_1.0.0 Biostrings_2.50.1
## [37] RJSONIO_1.3-1.1 rtracklayer_1.42.1
## [39] xfun_0.4 stringr_1.3.1
## [41] proto_1.0.0 ensembldb_2.6.3
## [43] zlibbioc_1.28.0 scales_1.0.0
## [45] BSgenome_1.50.0 VariantAnnotation_1.28.3
## [47] ProtGenerics_1.14.0 hms_0.4.2
## [49] SummarizedExperiment_1.12.0 RBGL_1.58.1
## [51] AnnotationFilter_1.6.0 yaml_2.2.0
## [53] curl_3.2 memoise_1.1.0
## [55] gridExtra_2.3 ggplot2_3.1.0
## [57] downloader_0.4 biomaRt_2.38.0
## [59] rpart_4.1-13 latticeExtra_0.6-28
## [61] stringi_1.2.4 checkmate_1.8.5
## [63] GenomicFeatures_1.34.1 BiocParallel_1.16.2
## [65] chron_2.3-53 rlang_0.3.0.1
## [67] pkgconfig_2.0.2 matrixStats_0.54.0
## [69] evaluate_0.12 lattice_0.20-38
## [71] purrr_0.2.5 bindr_0.1.1
## [73] GenomicAlignments_1.18.0 htmlwidgets_1.3
## [75] bit_1.1-14 tidyselect_0.2.5
## [77] plyr_1.8.4 magrittr_1.5
## [79] bookdown_0.8 R6_2.3.0
## [81] Hmisc_4.1-1 DelayedArray_0.8.0
## [83] DBI_1.0.0 gsubfn_0.7
## [85] pillar_1.3.0 foreign_0.8-71
## [87] survival_2.43-3 nnet_7.3-12
## [89] tibble_1.4.2 crayon_1.3.4
## [91] rmarkdown_1.10 progress_1.2.0
## [93] data.table_1.11.8 blob_1.1.1
## [95] digest_0.6.18 xtable_1.8-3
## [97] munsell_0.5.0 tcltk_3.5.1
The InterMine model could be accessed from the mine homepage by clicking the tab “QueryBuilder” and selecting the appropriate data type under “Select a Data Type to Begin a Query”:
Here we select Gene as the data type:
We could select Symbol and Chromosome->Primary Identifier by clicking Show on the right of them. Then click “Export XML” at the bottom right corner of the webpage:
The column names Gene.symbol and Gene.chromosome.primaryIdentifier are contained in the XML output: