View source: R/subAxt-methods.R
psubAxt | R Documentation |
Given two GRanges
objects, select the Axt
alignments
whose the target and query alignments are both within each pair of ranges.
psubAxt(x, targetSearch, querySearch)
x |
|
targetSearch,querySearch |
|
The ‘targetSearch’ and ‘querySearch’ have the coordinates relative to the positive strand. For each pair of the ranges, the alignments that lie within both the target and query range are kept.
A Axt
object.
Ge Tan
psubAxt
library(GenomicRanges) tAssemblyFn <- file.path(system.file("extdata", package="BSgenome.Drerio.UCSC.danRer10"), "single_sequences.2bit") qAssemblyFn <- file.path(system.file("extdata", package="BSgenome.Hsapiens.UCSC.hg38"), "single_sequences.2bit") axtFn <- file.path(system.file("extdata", package="CNEr"), "danRer10.hg38.net.axt") axt <- readAxt(axtFn, tAssemblyFn, qAssemblyFn) targetSearch <- GRanges(seqnames=c("chr6"), ranges=IRanges(start=c(24000000, 26900000), end=c(24060000, 26905000)), strand="+" ) querySearch <- GRanges(seqnames=c("chr7", "chr2"), ranges=IRanges(start=c(12577000, 241262700), end=c(12579000, 241268600)), strand="+" ) psubAxt(axt, targetSearch, querySearch)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.