View source: R/AnnotationGenome.R
createGeneAnnotation | R Documentation |
This function will create a gene annotation object that can be used for creating ArrowFiles or an ArchRProject, etc.
createGeneAnnotation(
genome = NULL,
TxDb = NULL,
OrgDb = NULL,
genes = NULL,
exons = NULL,
TSS = NULL,
annoStyle = NULL,
singleStrand = TRUE
)
genome |
A string that specifies the genome (ie "hg38", "hg19", "mm10", "mm9"). If |
TxDb |
A |
OrgDb |
An |
genes |
A |
exons |
A |
TSS |
A |
annoStyle |
annotation style to map between gene names and various gene identifiers e.g. "ENTREZID", "ENSEMBL". |
singleStrand |
A boolean for GenomicFeatures::genes( |
if (!require("TxDb.Hsapiens.UCSC.hg19.knownGene", quietly = TRUE)) BiocManager::install("TxDb.Hsapiens.UCSC.hg19.knownGene", update = FALSE)
if (!require("org.Hs.eg.db", quietly = TRUE)) BiocManager::install("org.Hs.eg.db", update = FALSE)
library(TxDb.Hsapiens.UCSC.hg19.knownGene)
library(org.Hs.eg.db)
# Get Txdb
TxDb <- TxDb.Hsapiens.UCSC.hg19.knownGene
# Get OrgDb
OrgDb <- org.Hs.eg.db
# Create Genome Annotation
geneAnno <- createGeneAnnotation(TxDb=TxDb, OrgDb=OrgDb)
# Also can create from a string if BSgenome exists
geneAnno <- createGeneAnnotation("hg19")
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.