getFirstStrandedFromGRL: Get first stranded GRanges feature per GRangesList

getFirstStrandedFromGRLR Documentation

Get first stranded GRanges feature per GRangesList

Description

Get first stranded GRanges feature per GRangesList

Usage

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.

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.

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.

See Also

Other jam ALE-specific RNA-seq functions: ale2violin(), tx2ale()

Other jam GRanges functions: addGRLgaps(), addGRgaps(), annotateGRLfromGRL(), annotateGRfromGR(), assignGRLexonNames(), closestExonToJunctions(), combineGRcoverage(), exoncov2polygon(), findOverlapsGRL(), flattenExonsBy(), getGRLgaps(), getGRcoverageFromBw(), getGRgaps(), grl2df(), jam_isDisjoint(), make_ref2compressed(), sortGRL(), spliceGR2junctionDF(), stackJunctions()

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)


jmw86069/splicejam documentation built on Nov. 4, 2024, 10:53 a.m.