reproduce_ranges: Reproduce the gene or exons used in the...

View source: R/reproduce_ranges.R

reproduce_rangesR Documentation

Reproduce the gene or exons used in the RangedSummarizedExperiment objects

Description

This function reproduces the gene or exon level information used for creating the RangedSummarizedExperiment-class objects provided by recount. The annotation is based on Gencode v25 with the gene-level information extracted with genes() (see transcripts with default arguments.

Usage

reproduce_ranges(level = "gene", db = "Gencode.v25")

Arguments

level

Either genes or exon. It specifies whether to return Gene or exon level information as a GRanges-class or GRangesList-class object respectively. The gene level information contains the width of the disjoint exons for that given gene which can be used to normalize the counts provided by recount. Can also be both in which case a 2 element list with the exon and the gene output is returned.

db

Either Gencode.v25 (default) or EnsDb.Hsapiens.v79. The default option reproduces the annotation used when creating recount. EnsDb.Hsapiens.v79 can be used for an alternative annotation as showcased in the recount vignette.

Details

For Gencode.v25, we use the comprehensive gene annotation (regions: CHR) from https://www.gencodegenes.org/releases/25.html (GRCh38.p7).

Note that gene symbols have changed over time. This answer in the Bioconductor support forum details how to obtain the latest gene symbol mappings: https://support.bioconductor.org/p/126148/#126173.

Value

Either a GRanges-class object like recount_genes or a GRangesList-class object like recount_exons.

Author(s)

Leonardo Collado-Torres

See Also

recount_genes, recount_exons, https://github.com/nellore, https://jhubiostatistics.shinyapps.io/recount/

Examples

## Not run: 
## Reproduce gene level information
genes <- reproduce_ranges()

## Compare against recount_genes
length(genes)
length(recount_genes)

## End(Not run)


leekgroup/recount documentation built on March 28, 2024, 10:48 a.m.