Title: | An R package to perform genomic liftover |
---|---|
Description: | The easylift package provides a convenient tool for genomic liftover operations between different genome assemblies. It seamlessly works with Bioconductor's GRanges objects and chain files from the UCSC Genome Browser, allowing for straightforward handling of genomic ranges across various genome versions. One noteworthy feature of easylift is its integration with the BiocFileCache package. This integration automates the management and caching of chain files necessary for liftover operations. Users no longer need to manually specify chain file paths in their function calls, reducing the complexity of the liftover process. |
Authors: | Abdullah Al Nahid [aut, cre] , Hervé Pagès [aut, rev], Michael Love [aut, rev] |
Maintainer: | Abdullah Al Nahid <[email protected]> |
License: | MIT + file LICENSE |
Version: | 1.5.0 |
Built: | 2024-10-31 18:10:13 UTC |
Source: | https://github.com/bioc/easylift |
This function takes a GRanges object with genomic coordinates in one genome assembly and lifts them to target genome assembly using a chain file.
easylift(x, to, chain, bfc)
easylift(x, to, chain, bfc)
x |
A GRanges object with genomic coordinates in the original assembly. |
to |
The target genome assembly (e.g., "hg38"). |
chain |
The path to the chain file containing the liftover mapping. Can be provided in gzipped or non-gzipped format. If omitted, the function will look in the default BiocFileCache for a properly named chain file. |
bfc |
A BiocFileCache object (optional), if not provided (most typically) the default location will be used. |
A GRanges object with lifted genomic coordinates.
liftOver
function from the rtracklayer
package,
which is the basis for easylift
.
# Lift over the coordinates of the first 10 genes in the hg19 assembly # to the hg38 assembly library(easylift) gr <- GRanges( seqname = Rle(c("chr1", "chr2"), c(100000, 100000)), ranges = IRanges(start = 1, end = 200000) ) # Here, "hg19" is the source genome genome(gr) <- "hg19" # Here, we use the `system.file()` function because the chain file is in the # package (however if you need to point to any other file on your machine, # just do 'chain <- "path/to/your/hg19ToHg38.over.chain.gz"'): chain <- system.file("extdata", "hg19ToHg38.over.chain.gz", package = "easylift") # Here, "hg38" is the target genome easylift(gr, "hg38", chain) # To use `BiocFileCache` for the chain file, add it to the cache as follows: chain_file <- "/path/to/your/hg19ToHg38.over.chain.gz" bfc <- BiocFileCache() # Add chain file to cache if already not available if (nrow(bfcquery(bfc, basename(chain_file))) == 0) bfcadd(bfc, chain_file) # Then, use it in `easylift` like this: easylift(gr, "hg38") # or gr |> easylift("hg38")
# Lift over the coordinates of the first 10 genes in the hg19 assembly # to the hg38 assembly library(easylift) gr <- GRanges( seqname = Rle(c("chr1", "chr2"), c(100000, 100000)), ranges = IRanges(start = 1, end = 200000) ) # Here, "hg19" is the source genome genome(gr) <- "hg19" # Here, we use the `system.file()` function because the chain file is in the # package (however if you need to point to any other file on your machine, # just do 'chain <- "path/to/your/hg19ToHg38.over.chain.gz"'): chain <- system.file("extdata", "hg19ToHg38.over.chain.gz", package = "easylift") # Here, "hg38" is the target genome easylift(gr, "hg38", chain) # To use `BiocFileCache` for the chain file, add it to the cache as follows: chain_file <- "/path/to/your/hg19ToHg38.over.chain.gz" bfc <- BiocFileCache() # Add chain file to cache if already not available if (nrow(bfcquery(bfc, basename(chain_file))) == 0) bfcadd(bfc, chain_file) # Then, use it in `easylift` like this: easylift(gr, "hg38") # or gr |> easylift("hg38")