gscores {GenomicScores} | R Documentation |
Functions to access genomic gscores through GScores
objects.
availableGScores() getGScores(x) ## S4 method for signature 'GScores,GenomicRanges' gscores(x, ranges, ...) ## S4 method for signature 'GScores,character' gscores(x, ranges, ...) ## S4 method for signature 'GScores' score(x, ..., simplify=TRUE)
x |
For |
ranges |
A |
... |
In the call to the
|
simplify |
Flag setting whether the result should be simplified to a
vector ( |
The function availableGScores()
shows genomic score sets available as
AnnotationHub
online resources.
The method gscores()
takes as first argument a GScores-class
object that can be loaded from an annotation package or from an
AnnotationHub
resource. These two possibilities are illustrated in the
examples below.
The arguments ref
and alt
serve two purposes. One, when there
are multiple scores per position, such as with CADD or M-CAP, and we want to
select a score matching a specific combination of reference and alternate
alleles. The other purpose is when the GScores
object x
is a
MafDb.*
package, then by providing ref
and alt
alelles
we will get separate frequencies for reference and alternate alleles. The
current lossy compression of these values yields a correct assignment for
biallelic variants in the corresponding MafDb.*
package and an
approximation for multiallelic ones.
The function availableGScores()
returns a character vector with the
names of the AnnotationHub
resources corresponding to different
available sets of genomic scores. The function getGScores()
return a
GScores
object. The method gscores()
returns a GRanges
object with the genomic scores in a metadata column called score
.
The method score()
returns a numeric vector with the genomic scores.
R. Castelo
Puigdevall, P. and Castelo, R. GenomicScores: seamless access to genomewide position-specific scores from R and Bioconductor. Bioinformatics, 18:3208-3210, 2018.
phastCons100way.UCSC.hg19
MafDb.1Kgenomes.phase1.hs37d5
## one genomic range of width 5 gr1 <- GRanges(seqnames="chr7", IRanges(start=117232380, width=5)) gr1 ## five genomic ranges of width 1 gr2 <- GRanges(seqnames="chr7", IRanges(start=117232380:117232384, width=1)) gr2 ## accessing genomic gscores from an annotation package if (require(phastCons100way.UCSC.hg19)) { library(GenomicRanges) gsco <- phastCons100way.UCSC.hg19 gsco gscores(gsco, gr1) score(gsco, gr1) gscores(gsco, gr2) populations(gsco) gscores(gsco, gr2, pop="DP2") } if (require(MafDb.1Kgenomes.phase1.hs37d5)) { mafdb <- MafDb.1Kgenomes.phase1.hs37d5 mafdb populations(mafdb) ## lookup allele frequencies for SNP rs1129038, located at 15:28356859, a ## SNP associated to blue and brown eye colors as reported by Eiberg et al. ## Blue eye color in humans may be caused by a perfectly associated founder ## mutation in a regulatory element located within the HERC2 gene ## inhibiting OCA2 expression. Human Genetics, 123(2):177-87, 2008 ## [http://www.ncbi.nlm.nih.gov/pubmed/18172690] gscores(mafdb, GRanges("15:28356859"), pop=populations(mafdb)) gscores(mafdb, "rs1129038", pop=populations(mafdb)) } ## accessing genomic scores from AnnotationHub resources ## Not run: availableGScores() gsco <- getGScores("phastCons100way.UCSC.hg19") gscores(gsco, gr1) ## End(Not run)