CTSS-class | R Documentation |
The CTSS
class represents CAGE transcription start sites (CTSS) at
single-nucleotide resolution, using GenomicRanges::UnstitchedGPos
as base
class. It is used by CAGEr for type safety.
The CTSS
constructor takes the same arguments as GenomicRanges::GPos
,
plus bsgenomeName
, and minus stitch
, which is hardcoded to FALSE
.
## S4 method for signature 'CTSS'
show(object)
## S4 method for signature 'CTSS'
initialize(.Object, ..., bsgenomeName = NULL)
CTSS(
seqnames = NULL,
pos = NULL,
strand = NULL,
...,
seqinfo = NULL,
seqlengths = NULL,
bsgenomeName = NULL
)
## S4 method for signature 'CTSS,GRanges'
coerce(from, to = "GRanges", strict = TRUE)
## S4 method for signature 'GRanges,CTSS'
coerce(from, to = "CTSS", strict = TRUE)
object |
See |
.Object |
See |
bsgenomeName |
String containing the name of a BSgenome package. |
seqnames , pos , strand , seqinfo , seqlengths , ... |
See the documentation
of |
from , to , strict |
See |
The genomeName
element of the metadata
slot is used to store the
name of the BSgenome package used when constructing the CAGEr
object.
Coercion from GRanges
to CTSS
loses information, but it seems
to be fine, since other coercions like as(1.2, "integer")
do the same.
Charles Plessy
# Convert an UnstitchedGPos object using the new() constructor.
gp <- GPos("chr1:2:-", stitch = FALSE)
ctss <- new("CTSS", gp, bsgenomeName = "BSgenome.Drerio.UCSC.danRer7")
genomeName(ctss)
# Create a new object using the CTSS() constructor.
CTSS("chr1", 2, "-", bsgenomeName = "BSgenome.Drerio.UCSC.danRer7")
# Coerce CTSS to GRanges
as(ctss, "GRanges")
# Coerce a GRanges object to CTSS using the as() method.
gr <- GRanges("chr1:1-10:-")
gr$seq <- "AAAAAAAAAA"
seqlengths(gr) <- 100
genome(gr) <- "foo"
as(gr, "CTSS")
identical(seqinfo(gr), seqinfo(as(gr, "CTSS")))
as(as(gr, "CTSS"), "CTSS") # Make sure it works twice in a row
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.