Axt-class: Class '"Axt"'

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

Description

The Axt S4 object to hold a axt file.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
## Constructors:
Axt(targetRanges=GRanges(), targetSeqs=DNAStringSet(),
    queryRanges=GRanges(), querySeqs=DNAStringSet(),
    score=integer(0), symCount=integer(0), names=NULL)

## Accessor-like methods:
## S4 method for signature 'Axt'
targetRanges(x)
## S4 method for signature 'Axt'
targetSeqs(x)
## S4 method for signature 'Axt'
queryRanges(x)
## S4 method for signature 'Axt'
querySeqs(x)
## S4 method for signature 'Axt'
score(x)
## S4 method for signature 'Axt'
symCount(x)
## ... and more (see Methods)

Arguments

targetRanges

Object of class "GRanges": The ranges of net alignments on reference genome.

targetSeqs

Object of class "DNAStringSet": The alignment sequences of reference genome.

queryRanges

Object of class "GRanges": The ranges of net alignments on query genome.

querySeqs

Object of class "DNAStringSet": The alignment sequences of query genome.

score

Object of class "integer": The alignment score.

symCount

Object of class "integer": The alignment length.

names

character(): the names of axt alignments.

x

Object of class "Axt": A Axt object.

Details

In ‘axt’ files and Axt object, the ‘targetRanges’ also have the alignments on positive strands. However, the ‘queryRanges’ can have alignments on negative strands, and the coordinates are based on negative strands, which is quite different from the convention in Bioconductor. To convert the coordinates of alignments on the negative strand to the positive strand, use normaliseStrand.

Methods

[

signature(x = "Axt", i = "ANY", j = "ANY"): Axt getter

c

signature(x = "Axt"): Axt concatenator.

length

signature(x = "Axt"): Get the number of alignments.

queryRanges

signature(x = "Axt"): Get the ranges of query genome.

querySeqs

signature(x = "Axt"): Get the alignment sequences of query genome.

score

signature(x = "Axt"): Get the alignment score.

symCount

signature(x = "Axt"): Get the alignment lengths.

targetRanges

signature(x = "Axt"): Get the ranges of reference genome.

targetSeqs

signature(x = "Axt"): Get the alignment sequences of reference genome.

Author(s)

Ge Tan

See Also

readAxt writeAxt subAxt fixCoordinates makeAxtTracks

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
  library(GenomicRanges)
  library(Biostrings)
  ## Constructor
  targetRanges <- GRanges(seqnames=c("chr1", "chr1", "chr2", "chr3"),
                          ranges=IRanges(start=c(1, 20, 2, 3),
                                         end=c(10, 25, 10, 10)),
                          strand="+")
  targetSeqs <- DNAStringSet(c("ATTTTATGTG", "GGGAAG", "GGGCTTTTG",
                               "TTGTGTAG"))
  queryRanges <- GRanges(seqnames=c("chr1", "chr10", "chr10", "chr20"),
                         ranges=IRanges(start=c(1, 25, 50, 5),
                                        end=c(10, 30, 58, 12)),
                         strand="+")
  querySeqs <- DNAStringSet(c("ATTTAAAGTG", "GGAAAA", "GGGCTCTGG",
                              "TTAAATAA"))
  score <- c(246L, 4422L, 5679L, 1743L)
  symCount <- c(10L, 6L, 9L, 8L)
  axt <- Axt(targetRanges=targetRanges, targetSeqs=targetSeqs, 
             queryRanges=queryRanges, querySeqs=querySeqs, 
             score=score, symCount=symCount)
  
  ## getters
  names(axt)
  length(axt)
  first(axt)
  last(axt)
  seqnames(axt)
  strand(axt)
  seqinfo(axt)
  
  ## Vector methods
  axt[1]
  
  ## List methods
  unlist(axt)
  
  ## Combining
  c(axt, axt)

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