conversions: Convert between GRanges and PyRanges objects

Bioc-PyRanges-conversionsR Documentation

Convert between GRanges and PyRanges objects

Description

Two utility functions for converting back and forth between GRanges and PyRanges objects.

Usage

makePyRangesFromGRanges(gr)

makeGRangesFromPyRanges(pyr)

Arguments

gr

A GRanges object, or, more generally, it can be any GenomicRanges derivative.

pyr

A PyRanges object.

Value

makePyRangesFromGRanges returns a PyRanges object.

makeGRangesFromPyRanges returns a GRanges object.

References

PyRanges: https://github.com/biocore-NTNU/pyranges

PyRanges' documentation: https://biocore-NTNU.github.io/pyranges/

See Also

  • GRanges objects defined in the GenomicRanges package.

  • What is the BiocPyRanges package about.

Examples

## ---------------------------------------------------------------------
## A SIMPLE EXAMPLE
## ---------------------------------------------------------------------

gr <- GRanges("chr1", IRanges(1:3, 20))
pyr <- makePyRangesFromGRanges(gr)
pyr  # PyRanges objects follow the 0-based start convention!!

makeGRangesFromPyRanges(pyr)

## ---------------------------------------------------------------------
## TYPICAL BiocPyRanges WORKFLOW
## ---------------------------------------------------------------------

gr1 <- GRanges("chr1", IRanges(c(1:3, 15, 35, 51), c(20:22, 30, 50, 60)),
               strand=c("+", "+", "-", "+", "+", "+"))
gr1

## The typical BiocPyRanges workflow has 3 steps:

## STEP 1 - Go to PyRanges space (i.e. turn your GRange object into a
##          PyRanges object):
pyr1 <- makePyRangesFromGRanges(gr1)
pyr1

## STEP 2 - Operate in PyRanges space (i.e. perform some operation
##          provided by the PyRanges module on your PyRanges object):
pyr2 <- pyr1$merge()  # equivalent of GenomicRanges::reduce()

## STEP 3 - Come back to R space (i.e. turn result back into GRanges
##          object):
gr2 <- makeGRangesFromPyRanges(pyr2)

## Sanity check:
stopifnot(identical(reduce(gr1, ignore.strand=TRUE), gr2))

## ---------------------------------------------------------------------
## PERFORMANCE
## ---------------------------------------------------------------------

## PyRanges is very fast!
random_GRanges <- function(n, width) {
    gr_starts <- sample(50e6, n, replace=TRUE)
    chroms <- c("chr1", "chr2", "chr3")
    gr_chrom <- sample(factor(chroms, levels=chroms), n, replace=TRUE)
    gr_strand <- sample(levels(strand()), n, replace=TRUE)
    GRanges(gr_chrom, IRanges(gr_starts, width=width, names=seq_len(n)),
            strand=gr_strand)
}

set.seed(123)
gr1 <- random_GRanges(3e6, width=40)  # 3 million random ranges

system.time(gr2a <- reduce(gr1, ignore.strand=TRUE))

pyr1 <- makePyRangesFromGRanges(gr1)
system.time(pyr2 <- pyr1$merge())
gr2b <- makeGRangesFromPyRanges(pyr2)

## Sanity check:
stopifnot(identical(gr2a, gr2b))

hpages/BiocPyRanges documentation built on Aug. 21, 2023, 10:11 a.m.