Skip to contents

Modify one or more sequences to include Insertions or Deletions

Usage

indelcator(x, indels, ...)

# S4 method for class 'XString,GRanges'
indelcator(x, indels, exons, alt_col = "ALT", ...)

# S4 method for class 'DNAStringSet,GRanges'
indelcator(x, indels, alt_col = "ALT", mc.cores = 1, verbose = TRUE, ...)

# S4 method for class 'BSgenome,GRanges'
indelcator(x, indels, alt_col = "ALT", mc.cores = 1, names, ...)

Arguments

x

Sequence of class XString

indels

GRanges object with InDel locations and the alternate allele

...

Passed to parallel::mclapply

exons

GRanges object containing exon structure for x

alt_col

Column containing the alternate allele

mc.cores

Number of cores to use when calling parallel::mclapply internally

verbose

logical(1) Print all messages

names

passed to BSgenome::getSeq when x is a BSgenome object

Value

A DNAStringSet or XString object (See Details)

Details

This is a lower-level function relied on by both transmogrify() and genomogrify().

Takes an Biostrings::XString or Biostrings::XStringSet object and modifies the sequence to incorporate InDels. The expected types of data determine the behaviour, with the following expectations describing how the function will incorporate data

Input Data TypeExons RequiredUse CaseReturned
XStringYModify a Reference TranscriptomeXString
DNAStringSetNModify a Reference GenomeDNAStringSet
BSgenomeNModify a Reference GenomeDNAStringSet

Examples

## Start with a DNAStringSet
library(GenomicRanges)
seq <- DNAStringSet(c(seq1 = "AATCTGCGC"))
## Define an Insertion
var <- GRanges("seq1:1")
var$ALT <- "AAA"
seq
#> DNAStringSet object of length 1:
#>     width seq                                               names               
#> [1]     9 AATCTGCGC                                         seq1
indelcator(seq, var)
#> Updating seq1; Original length: 9
#> ; Updated length: 11
#> DNAStringSet object of length 1:
#>     width seq                                               names               
#> [1]    11 AAAATCTGCGC                                       seq1

## To modify a single transcript
library(GenomicFeatures)
#> Loading required package: AnnotationDbi
#> Loading required package: Biobase
#> Welcome to Bioconductor
#> 
#>     Vignettes contain introductory material; view with
#>     'browseVignettes()'. To cite Bioconductor, see
#>     'citation("Biobase")', and for packages 'citation("pkgname")'.
ex <- GRanges(c("seq1:1-3:+", "seq1:7-9:+"))
orig <- extractTranscriptSeqs(seq, GRangesList(tx1 = ex))[["tx1"]]
orig
#> 6-letter DNAString object
#> seq: AATCGC
indelcator(orig, var, exons = ex)
#> 8-letter DNAString object
#> seq: AAAATCGC