ranges-anchor: Anchored Ranges objects

anchorR Documentation

Anchored Ranges objects

Description

The GRangesAnchored class and the IRangesAnchored class allow components of a GRanges or IRanges (start, end, center) to be held fixed.

Usage

anchor(x)

unanchor(x)

anchor_start(x)

anchor_end(x)

anchor_center(x)

anchor_centre(x)

anchor_3p(x)

anchor_5p(x)

Arguments

x

a Ranges object

Details

Anchoring will fix a Ranges start, end, or center positions, so these positions will remain the same when performing arithimetic. For GRanges objects, the function (anchor_3p()) will fix the start for the negative strand, while anchor_5p() will fix the end for the positive strand. Anchoring modifies how arithmetic is performed, for example modifying the width of a range with set_width() or stretching a range with stretch(). To remove anchoring use unanchor().

Value

a RangesAnchored object which has the same appearance as a regular Ranges object but with an additional slot displaying an anchor.

Constructors

Depending on how you want to fix the components of a Ranges, there are five ways to construct a RangesAnchored class. Here x is either an IRanges or GRanges object.

  • anchor_start(x)Fix the start coordinates

  • anchor_end(x)Fix the end coordinates

  • anchor_center(x)Fix the center coordinates

  • anchor_3p(x)On the negative strand fix the start coordinates, and for positive or unstranded ranges fix the end coordinates.

  • anchor_5p(x)On the positive or unstranded ranges fix the start coordinates, coordinates and for negative stranded ranges fix the end coordinates.

Accessors

To see what has been anchored use the function anchor. This will return a character vector containing a valid anchor. It will be set to one of c("start", "end", "center") for an IRanges object or one of c("start", "end", "center", "3p", "5p") for a GRanges object.

See Also

mutate, stretch

Examples

df <- data.frame(start = 1:10, width = 5)
rng <- as_iranges(df)
rng_by_start <- anchor_start(rng)
rng_by_start
anchor(rng_by_start)
mutate(rng_by_start, width = 3L)
grng <- as_granges(df,
                   seqnames = "chr1",
                   strand = c(rep("-", 5), rep("+", 5)))
rng_by_5p <- anchor_5p(grng)
rng_by_5p
mutate(rng_by_5p, width = 3L)


sa-lee/plyranges documentation built on April 15, 2024, 12:25 p.m.