snps.from.rsid: Import SNPs from rsid for use in motifbreakR

View source: R/locus.from.rsid.R

snps.from.rsidR Documentation

Import SNPs from rsid for use in motifbreakR

Description

Import SNPs from rsid for use in motifbreakR

Usage

snps.from.rsid(
  rsid = NULL,
  dbSNP = NULL,
  search.genome = NULL,
  biomart.dataset = NULL
)

Arguments

rsid

Character; a character vector of rsid values from dbSNP

dbSNP

an object of class SNPlocs to lookup rsids; see availible.SNPs in injectSNPs to check for availible SNPlocs

search.genome

an object of class BSgenome for the species you are interrogating; see available.genomes for a list of species.

biomart.dataset

a Mart object from useEnsembl specifying the snps biomart, which dataset i.e., hsapiens_snp, and which version i.e., 111 or GRCm39. This will override SNPlocs and must be compatible with search.genome selection, and will query from biomaRt, which may be considerably faster than a lookup from a SNPlocs object.

Details

snps.from.rsid take an rsid, or character vector of rsids and generates the required object to input into motifbreakR

Value

a GRanges object containing:

SNP_id

The rsid of the snp with the "rs" portion stripped

alleles_as_ambig

THE IUPAC ambiguity code between the reference and alternate allele for this SNP

REF

The reference allele for the SNP

ALT

The alternate allele for the SNP

See Also

See motifbreakR for analysis; See snps.from.file for an alternate method for generating a list of variants.

Examples

 library(BSgenome.Hsapiens.UCSC.hg19)
 library(SNPlocs.Hsapiens.dbSNP155.GRCh37)
 snps.file <- system.file("extdata", "pca.enhancer.snps", package = "motifbreakR")
 snps <- as.character(read.table(snps.file)[,1])
 snps.mb <- snps.from.rsid(snps[1],
                           dbSNP = SNPlocs.Hsapiens.dbSNP155.GRCh37,
                           search.genome = BSgenome.Hsapiens.UCSC.hg19)

 ## alternatively using biomaRt

 library(biomaRt)
 library(BSgenome.Hsapiens.UCSC.hg38)
 ensembl_snp <- useEnsembl(biomart = "snps",
                           dataset = "hsapiens_snp",
                           version = "112")
 snps.mb <- snps.from.rsid(snps,
                           biomart.dataset = ensembl_snp,
                           search.genome = BSgenome.Hsapiens.UCSC.hg38)


Simon-Coetzee/motifBreakR documentation built on Aug. 6, 2024, 5:17 a.m.