Description Usage Arguments Details Value References See Also Examples
This method retrieves and stores GO annotations for the organism of interest from one of genomic ressource database (Bioconductor, EntrezGene, Ensembl, Uniprot).
1 2 3 4 |
id |
identifiant corresponding to the organism of interest.
This id name is referenced in the first column of the database
used (see |
object |
a required |
ortholog |
|
This method uses a genomic_ressource-class
object to retrieve
GO annotations for the organism of interest.
The stored annotations are structured in 3 slots corresponding to the 3 GO categories: MF (Molecular Function),
BP (Biological Process), and CC (Cellular Component). Each slot contains GO terms with
associated evidence code.
The genomic_ressource-class
object is created by one of the four available methods:
Bioconductor2GO
, EntrezGene2GO
,
Ensembl2GO
, or Uniprot2GO
.
In the case of vertebrates, setting ortholog
argument to TRUE
is required if you need to add GO terms with experimental
evidence codes from orthologs genes
when using EntrezGene2GO
method. To display organisms supported by NCBI EntrezGene orthologs pipeline,
set the arguments id=NULL
and ortholog=TRUE
.
This approch is highly similar to the strategy developed by Uniprot-GOA consortium for the Electronic Annotation Method using
Ensembl Compara.
annotate
produces an object of gene2GO-class
required by build_GO_SS
method.
Durinck S, Spellman P, Birney E and Huber W (2009). Mapping identifiers for the integration of genomic datasets with the R/Bioconductor package biomaRt. Nature Protocols, 4, pp. 1184-1191.
Durinck S, Moreau Y, Kasprzyk A, Davis S, De Moor B, Brazma A and Huber W (2005). BioMart and Bioconductor: a powerful link between biological databases and microarray data analysis. Bioinformatics, 21, pp. 3439-3440.
Fong, JH, Murphy, TD, Pruitt, KD (2013). Comparison of RefSeq protein-coding regions in human and vertebrate genomes. BMC Genomics, 14:654.
Henrik Bengtsson (2016). R.utils: Various Programming Utilities. R package version 2.5.0. https://CRAN.R-project.org/package=R.utils.
Herve Pages, Marc Carlson, Seth Falcon and Nianhua Li (2017). AnnotationDbi: Annotation Database Interface. R package version 1.38.0.
Matt Dowle and Arun Srinivasan (2017). data.table: Extension of data.frame. R package version 1.10.4. https://CRAN.R-project.org/package=data.table.
Other genomic_ressource:
Bioconductor2GO()
,
Custom2GO()
,
Ensembl2GO()
,
EntrezGene2GO()
,
Uniprot2GO()
,
available_organisms()
,
genomic_ressource-class
,
taxonomy()
Other GO_terms:
GOcount()
,
GOterms_heatmap()
,
create_topGOdata()
,
gene2GO-class
,
merge_enrich_terms()
,
runfgsea()
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 | ## Not run:
## load Mus musculus (mouse) GO annotations
# from Bioconductor
Bioconductor<-ViSEAGO::Bioconductor2GO()
myGENE2GO<-ViSEAGO::annotate(
id="org.Mm.eg.db",
object=Bioconductor
)
# from EntrezGene
EntrezGene<-ViSEAGO::EntrezGene2GO()
myGENE2GO<-ViSEAGO::annotate(
id="10090",
object=EntrezGene
)
# from EntrezGene
Ensembl<-ViSEAGO::Ensembl2GO()
myGENE2GO<-ViSEAGO::annotate(
id="mmusculus_gene_ensembl",
object=Ensembl
)
# from Uniprot
Uniprot<-ViSEAGO::Uniprot2GO()
myGENE2GO<-ViSEAGO::annotate(
id="mouse",
object=Uniprot
)
## from Custom GO annotation file
Custom<-ViSEAGO::Custom2GO(system.file("extdata/customfile.txt",package = "ViSEAGO"))
myGENE2GO<-ViSEAGO::annotate(
id="myspecies1",
object=Custom
)
## specific options for EntrezGene database
# Chicken GO annotations without adding orthologs
EntrezGene<-ViSEAGO::EntrezGene2GO()
myGENE2GO<-ViSEAGO::annotate(
id="9031",
object=EntrezGene
)
# Chicken GO annotation with the add of orthologs GO annotations
EntrezGene<-ViSEAGO::EntrezGene2GO()
myGENE2GO<-ViSEAGO::annotate(
id="9031",
object=EntrezGene,
ortholog=TRUE
)
# display organisms supported by NCBI EntrezGene orthologs pipeline
EntrezGene<-ViSEAGO::EntrezGene2GO()
ViSEAGO::annotate(
id="NULL",
object=EntrezGene,
ortholog=TRUE
)
## End(Not run)
|
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.