annotateCTSS: Annotate and compute summary statistics

annotateCTSSR Documentation

Annotate and compute summary statistics

Description

annotateCTSS annotates the CTSS of a CAGEexp object and computes annotation statistics.

annotateConsensusClusters annotates the consensus clusters of a CAGEr object.

Usage

annotateCTSS(object, annot, upstream = 500, downstream = 500)

## S4 method for signature 'CAGEexp,GRanges'
annotateCTSS(object, annot, upstream = 500, downstream = 500)

## S4 method for signature 'CAGEexp,TxDb'
annotateCTSS(object, annot)

annotateTagClusters(object, annot, upstream = 500, downstream = 500)

## S4 method for signature 'CAGEexp,GRanges'
annotateTagClusters(object, annot, upstream = 500, downstream = 500)

## S4 method for signature 'CAGEexp,TxDb'
annotateTagClusters(object, annot)

annotateConsensusClusters(object, annot, upstream = 500, downstream = 500)

## S4 method for signature 'CAGEexp,GRanges'
annotateConsensusClusters(object, annot, upstream = 500, downstream = 500)

## S4 method for signature 'CAGEexp,TxDb'
annotateConsensusClusters(object, annot)

Arguments

object

CAGEexp object.

annot

A GRanges or a TxDb object representing the genome annotation. See details for the GRanges object.

upstream

Number of bases upstream the start of the transcript models to be considered as part of the promoter region.

downstream

Number of bases downstream the start of the transcript models to be considered as part of the promoter region.

Details

If the annotation is a GRanges, gene names will be extracted from the gene_name metadata, the transcript_type metadata will be used to filter out entries that do not have promoters (such as immunogloblulin VDJ segments), and the type metadata is used to extract positions of introns and exons.

Value

annotateCTSS returns the input object with the following modifications:

  • The Genomic Ranges of the tagCountMatrix experiment gains an annotation metadata column, with levels such as promoter, exon, intron and unknown. If the annotation has a gene_name metadata, then a genes column is also added, with gene symbols from the annotation.

  • The sample metadata gets new columns, indicating total counts in each of the annotation levels. If the annotation has a gene_name metadata, then a genes column is added to indicate the number of different gene symbols detected.

annotateTagClusters returns the input object with the same modifications as above.

annotateConsensusClusters returns the input object with the same modifications as above.

Author(s)

Charles Plessy

See Also

CTSStoGenes, and the exampleZv9_annot example data.

Other CAGEr object modifiers: CTSStoGenes(), CustomConsensusClusters(), aggregateTagClusters(), cumulativeCTSSdistribution(), distclu(), getCTSS(), normalizeTagCount(), paraclu(), quantilePositions(), quickEnhancers(), resetCAGEexp(), summariseChrExpr()

Other CAGEr annotation functions: plotAnnot(), ranges2annot(), ranges2genes(), ranges2names()

Examples

annotateCTSS(exampleCAGEexp, exampleZv9_annot)
colData(exampleCAGEexp)

exampleCAGEexp <- annotateTagClusters(exampleCAGEexp, exampleZv9_annot)
tagClustersGR(exampleCAGEexp, 1)

annotateConsensusClusters(exampleCAGEexp, exampleZv9_annot)
consensusClustersGR(exampleCAGEexp)


charles-plessy/CAGEr documentation built on Oct. 27, 2024, 10:11 p.m.