subAxt-methods: Subset an 'Axt' object

Description Usage Arguments Details Value Author(s) See Also Examples

Description

A ‘subAxt’ method for extracting a set of alignments from an Axt object.

Usage

1
  subAxt(x, chr, start, end, select=c("target", "query"), qSize=NULL)

Arguments

x

An object of Axt.

chr

An object of character containing the names of the sequences in 'x' where to get the alignments from, or a GRanges object where 'start' and 'end' are missing. In the case of GRanges, the strand information is ignored.

start, end

An object of integer() or missing. These ranges should be based on the positive strand. When select is "query", the reverse complement alignments which lay inside this range will also be selected.

select

When select is ‘target’, the subset criteria are applied on target alignments in Axt. When select is ‘query’, the subset criteria are applied on query alignments in Axt.

qSize

integer(n): When select is ‘query’, ‘qSize’ must exist in ‘x’ or can be provided as a vector of chromosome lengths.

Details

Usually when we want to subset some axts from a Axt object, we care about all the axts within a certain range. The axts can come from the axt file with chr as reference (i.e., target sequence), or the axt file with chr as query sequence. When the chr is query sequence, it can be on the negative strand. Hence, the size of chromosome is necessary to convert the search range to a range on negative strand coordinate.

When one Axt alignment partially overlaps the range, the whole Axt alignment will be extracted.

Value

An extracted Axt object is returned.

Author(s)

Ge Tan

See Also

psubAxt

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
  library(GenomicRanges)
  library(rtracklayer)
  
  ## Prepare the axt object
  tAssemblyFn <- file.path(system.file("extdata",
                           package="BSgenome.Hsapiens.UCSC.hg38"),
                           "single_sequences.2bit")
  qAssemblyFn <- file.path(system.file("extdata",
                           package="BSgenome.Drerio.UCSC.danRer10"),
                           "single_sequences.2bit")
  axtFilesHg38DanRer10 <- file.path(system.file("extdata", package="CNEr"), 
                                  "hg38.danRer10.net.axt")
  axtHg38DanRer10 <- readAxt(axtFilesHg38DanRer10, tAssemblyFn, qAssemblyFn)
  
  ## "character", "integer", "integer" on "target" sequence
  subAxt(axtHg38DanRer10, chr="chr1", start=148165963L, end=222131835L,
         select="target")
  
  ## "GRanges"" on "target" sequence
  searchGRanges <- GRanges(seqnames="chr1",
                           ranges=IRanges(start=148165963L,
                                          end=222131835L),
                           strand="+")
  subAxt(axtHg38DanRer10, searchGRanges, select="target")
  
  ## multiple "character", "integer", "integer" on "target" sequence
  subAxt(axtHg38DanRer10, chr=c("chr1", "chr13"), 
         start=c(148165963L, 94750629L),
         end=c(222131835L, 94966991L), select="target")
  
  ## "character" only on "target" sequence
  subAxt(axtHg38DanRer10, chr="chr1", select="target")
  
  ## GRanges on "query" sequence
  searchGRanges <- GRanges(seqnames="chr6",
                           ranges=IRanges(start=25825774,
                                          end=26745499),
                           strand="+")
  subAxt(axtHg38DanRer10, searchGRanges, select="query")

CNEr documentation built on Nov. 8, 2020, 5:36 p.m.