CTSS-class: CAGE Transcription Start Sites

CTSS-classR Documentation

CAGE Transcription Start Sites

Description

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.

The CTSS.chr class represents a CTSS object that is guaranteed to be only on a single chromosome. It is used internally by CAGEr for type-safe polymorphic dispatch.

Usage

## 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)

Arguments

object

See methods::show

.Object

See methods::new

bsgenomeName

String containing the name of a BSgenome package.

seqnames, pos, strand, seqinfo, seqlengths, ...

See the documentation of GenomicRanges::GPos for further details.

from, to, strict

See methods::coerce.

Details

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.

Author(s)

Charles Plessy

Examples

# 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

# (internal use) Transform CTSS to CTSS.chr object
ctss.chr <- as(CTSScoordinatesGR(exampleCAGEexp), "CTSS.chr")

charles-plessy/CAGEr documentation built on Nov. 4, 2023, 11:57 a.m.