Get first stranded GRanges feature per GRangesList

getFirstStrandedFromGRL(
  grl,
  method = c("direct", "endoapply", "flank"),
  verbose = FALSE,
  ...
)

Arguments

grl

GRangesList

method

character value in c("direct", "endoapply", "flank") representing which method to use to define the first feature. The "direct" method uses IRanges::heads() or IRanges::tails() depending upon the strandedness; the "endoapply" method uses S4Vectors::endoapply() to iterate each GRangesList element, then returns the result from head() or tail(). The "flank" method is intended to be equivalent but uses IRanges::heads() or IRanges::tails() then GenomicRanges::flank(), which is vectorized. The "flank" method is deprecated and "direct" is recommended as the preferred vectorized replacement.

verbose

logical indicating whether to print verbose output.

...

additional arguments are ignored.

Value

GRangesList containing only the first stranded GRanges feature per input GRangesList element. When method="flank" the output contains no metadata values, but this method is deprecated.

Details

This function returns the first feature per GRangesList, relative to the strand of the GRangesList. It assumes each GRangesList element has only one strand and one seqname, and will stop otherwise.

See also

Examples

gr12 <- GenomicRanges::GRanges(seqnames=rep("chr1", 9), ranges=IRanges::IRanges( start=c(100, 200, 400, 500, 300, 100, 200, 400, 600), width=c(50,50,50, 50,50,50, 50,50,50) ), strand=rep(c("+", "-", "+"), c(3,3,3)), gene_name=rep(c("GeneA", "GeneB", "GeneC"), each=3) ) grl <- GenomicRanges::split(gr12, GenomicRanges::values(gr12)$gene_name) getFirstStrandedFromGRL(grl)
#> Warning: failed to set names on the unlisted CompressedRleList object
#> GRangesList object of length 3: #> $GeneA #> GRanges object with 1 range and 1 metadata column: #> seqnames ranges strand | gene_name #> <Rle> <IRanges> <Rle> | <character> #> [1] chr1 100-149 + | GeneA #> ------- #> seqinfo: 1 sequence from an unspecified genome; no seqlengths #> #> $GeneB #> GRanges object with 1 range and 1 metadata column: #> seqnames ranges strand | gene_name #> <Rle> <IRanges> <Rle> | <character> #> [1] chr1 500-549 - | GeneB #> ------- #> seqinfo: 1 sequence from an unspecified genome; no seqlengths #> #> $GeneC #> GRanges object with 1 range and 1 metadata column: #> seqnames ranges strand | gene_name #> <Rle> <IRanges> <Rle> | <character> #> [1] chr1 200-249 + | GeneC #> ------- #> seqinfo: 1 sequence from an unspecified genome; no seqlengths #>