| Title: | Disk-backed Representation of ONT bedMethyl Files |
|---|---|
| Description: | Bioconductor-native infrastructure for handling large nanoporetech modkit bedMethyl pileup files from ONT data using HDF5Array and DelayedArray. |
| Authors: | Vasileios Lemonidis [aut, cre, cph] (ORCID: <https://orcid.org/0000-0002-6446-3536>), Center for Oncological Research, University of Antwerp [cph, fnd], Stichting Tegen Kanker [fnd] |
| Maintainer: | Vasileios Lemonidis <[email protected]> |
| License: | GPL (>= 2) |
| Version: | 1.1.0 |
| Built: | 2026-05-30 09:37:09 UTC |
| Source: | https://github.com/bioc/RBedMethyl |
Subset an RBedMethyl object by integer, logical, or GRanges index.
## S4 method for signature 'RBedMethyl,missing,missing,missing' x[i, j, ..., drop = TRUE]## S4 method for signature 'RBedMethyl,missing,missing,missing' x[i, j, ..., drop = TRUE]
x |
An |
i |
Integer, logical, or |
j |
Unused. |
... |
Unused. |
drop |
Unused. |
A filtered RBedMethyl object.
Returns a data.frame describing retrievable bedMethyl fields and their types.
bedMethylFields()bedMethylFields()
A data.frame with columns field, type, and description.
bedMethylFields()bedMethylFields()
Compute per-site methylation fraction for an RBedMethyl object.
Requires the mod_reads assay to be loaded.
## S4 method for signature 'RBedMethyl,missing' beta(a, b)## S4 method for signature 'RBedMethyl,missing' beta(a, b)
a |
An |
b |
Unused, kept for |
Numeric vector of per-site methylation fractions.
Filter an RBedMethyl object by minimum coverage.
filterByCoverage(x, min_cov)filterByCoverage(x, min_cov)
x |
An |
min_cov |
Minimum coverage threshold. |
A filtered RBedMethyl object.
lines <- c( paste("chr1", 0, 1, "m", 0, "+", 0, 1, 0, 10, 0.5, 5, 5, 0, 0, 0, 0, 0, sep = "\t"), paste("chr1", 10, 11, "m", 0, "+", 10, 11, 0, 20, 0.25, 5, 15, 0, 0, 0, 0, 0, sep = "\t") ) tmp <- tempfile(fileext = ".bed") writeLines(lines, tmp) bm <- readBedMethyl(tmp, mod = "m", fields = c("coverage", "pct", "mod_reads")) bm2 <- filterByCoverage(bm, min_cov = 15) length(RBedMethyl::beta(bm2))lines <- c( paste("chr1", 0, 1, "m", 0, "+", 0, 1, 0, 10, 0.5, 5, 5, 0, 0, 0, 0, 0, sep = "\t"), paste("chr1", 10, 11, "m", 0, "+", 10, 11, 0, 20, 0.25, 5, 15, 0, 0, 0, 0, 0, sep = "\t") ) tmp <- tempfile(fileext = ".bed") writeLines(lines, tmp) bm <- readBedMethyl(tmp, mod = "m", fields = c("coverage", "pct", "mod_reads")) bm2 <- filterByCoverage(bm, min_cov = 15) length(RBedMethyl::beta(bm2))
Disk-backed representation of nanoporetech modkit bedMethyl data from ONT sequencing.
assaysA SimpleList of assay arrays.
chrom_levelsCharacter vector of chromosome names.
strand_levelsCharacter vector of strand levels.
chr_indexMatrix of chromosome row ranges (start/end).
indexInteger vector of active row indices.
modModification code.
Create an RBedMethyl object backed by HDF5Array from a nanoporetech
modkit bedMethyl file (headerless).
readBedMethyl( bedmethyl, mod = "m", chunk_size = 5e+06, h5file = NULL, check_sorted = TRUE, fields = c("coverage", "mod_reads") )readBedMethyl( bedmethyl, mod = "m", chunk_size = 5e+06, h5file = NULL, check_sorted = TRUE, fields = c("coverage", "mod_reads") )
bedmethyl |
Path to a nanoporetech modkit bedMethyl file (optionally gzipped). |
mod |
Modification code to retain ( |
chunk_size |
Reserved for future use. |
h5file |
Path to the HDF5 file to create. Defaults to a deterministic
path in |
check_sorted |
Logical, check that records are sorted by chrom and chromStart. |
fields |
Character vector of numeric fields to load. Defaults to |
An RBedMethyl object.
lines <- c( paste("chr1", 0, 1, "m", 0, "+", 0, 1, 0, 10, 0.5, 5, 5, 0, 0, 0, 0, 0, sep = "\t"), paste("chr1", 10, 11, "m", 0, "+", 10, 11, 0, 20, 0.25, 5, 15, 0, 0, 0, 0, 0, sep = "\t") ) tmp <- tempfile(fileext = ".bed") writeLines(lines, tmp) bm <- readBedMethyl(tmp, mod = "m", fields = c("coverage", "pct", "mod_reads")) bmlines <- c( paste("chr1", 0, 1, "m", 0, "+", 0, 1, 0, 10, 0.5, 5, 5, 0, 0, 0, 0, 0, sep = "\t"), paste("chr1", 10, 11, "m", 0, "+", 10, 11, 0, 20, 0.25, 5, 15, 0, 0, 0, 0, 0, sep = "\t") ) tmp <- tempfile(fileext = ".bed") writeLines(lines, tmp) bm <- readBedMethyl(tmp, mod = "m", fields = c("coverage", "pct", "mod_reads")) bm
Subset an RBedMethyl object using a predicate over an assay.
subsetBy(x, column, FUN)subsetBy(x, column, FUN)
x |
An |
column |
Assay name to filter on (must be loaded). |
FUN |
Predicate function returning a logical vector. |
A filtered RBedMethyl object.
lines <- c( paste("chr1", 0, 1, "m", 0, "+", 0, 1, 0, 10, 0.5, 5, 5, 0, 0, 0, 0, 0, sep = "\t"), paste("chr1", 10, 11, "m", 0, "+", 10, 11, 0, 20, 0.25, 5, 15, 0, 0, 0, 0, 0, sep = "\t") ) tmp <- tempfile(fileext = ".bed") writeLines(lines, tmp) bm <- readBedMethyl(tmp, mod = "m", fields = c("coverage", "pct", "mod_reads")) bm2 <- subsetBy(bm, "coverage", function(v) v >= 15) length(RBedMethyl::beta(bm2))lines <- c( paste("chr1", 0, 1, "m", 0, "+", 0, 1, 0, 10, 0.5, 5, 5, 0, 0, 0, 0, 0, sep = "\t"), paste("chr1", 10, 11, "m", 0, "+", 10, 11, 0, 20, 0.25, 5, 15, 0, 0, 0, 0, 0, sep = "\t") ) tmp <- tempfile(fileext = ".bed") writeLines(lines, tmp) bm <- readBedMethyl(tmp, mod = "m", fields = c("coverage", "pct", "mod_reads")) bm2 <- subsetBy(bm, "coverage", function(v) v >= 15) length(RBedMethyl::beta(bm2))
Subset an RBedMethyl object by one or more chromosomes.
subsetByChromosomes(x, chr)subsetByChromosomes(x, chr)
x |
An |
chr |
Character vector of chromosome names. |
A filtered RBedMethyl object.
lines <- c( paste("chr1", 0, 1, "m", 0, "+", 0, 1, 0, 10, 0.5, 5, 5, 0, 0, 0, 0, 0, sep = "\t"), paste("chr2", 10, 11, "m", 0, "+", 10, 11, 0, 20, 0.25, 5, 15, 0, 0, 0, 0, 0, sep = "\t") ) tmp <- tempfile(fileext = ".bed") writeLines(lines, tmp) bm <- readBedMethyl(tmp, mod = "m", fields = c("coverage", "pct", "mod_reads")) bm2 <- subsetByChromosomes(bm, c("chr1")) length(RBedMethyl::beta(bm2))lines <- c( paste("chr1", 0, 1, "m", 0, "+", 0, 1, 0, 10, 0.5, 5, 5, 0, 0, 0, 0, 0, sep = "\t"), paste("chr2", 10, 11, "m", 0, "+", 10, 11, 0, 20, 0.25, 5, 15, 0, 0, 0, 0, 0, sep = "\t") ) tmp <- tempfile(fileext = ".bed") writeLines(lines, tmp) bm <- readBedMethyl(tmp, mod = "m", fields = c("coverage", "pct", "mod_reads")) bm2 <- subsetByChromosomes(bm, c("chr1")) length(RBedMethyl::beta(bm2))
Subset an RBedMethyl object by genomic interval.
subsetByRegion(x, chr, start, end)subsetByRegion(x, chr, start, end)
x |
An |
chr |
Chromosome name. |
start |
Region start (0-based, half-open). |
end |
Region end. |
A filtered RBedMethyl object.
lines <- c( paste("chr1", 0, 1, "m", 0, "+", 0, 1, 0, 10, 0.5, 5, 5, 0, 0, 0, 0, 0, sep = "\t"), paste("chr1", 10, 11, "m", 0, "+", 10, 11, 0, 20, 0.25, 5, 15, 0, 0, 0, 0, 0, sep = "\t") ) tmp <- tempfile(fileext = ".bed") writeLines(lines, tmp) bm <- readBedMethyl(tmp, mod = "m", fields = c("coverage", "pct", "mod_reads")) bm2 <- subsetByRegion(bm, "chr1", 0, 5) length(RBedMethyl::beta(bm2))lines <- c( paste("chr1", 0, 1, "m", 0, "+", 0, 1, 0, 10, 0.5, 5, 5, 0, 0, 0, 0, 0, sep = "\t"), paste("chr1", 10, 11, "m", 0, "+", 10, 11, 0, 20, 0.25, 5, 15, 0, 0, 0, 0, 0, sep = "\t") ) tmp <- tempfile(fileext = ".bed") writeLines(lines, tmp) bm <- readBedMethyl(tmp, mod = "m", fields = c("coverage", "pct", "mod_reads")) bm2 <- subsetByRegion(bm, "chr1", 0, 5) length(RBedMethyl::beta(bm2))
Subset an RBedMethyl object by overlaps with a GRanges.
## S4 method for signature 'RBedMethyl,GRanges,missing,missing' subsetByRegion(x, chr, start, end)## S4 method for signature 'RBedMethyl,GRanges,missing,missing' subsetByRegion(x, chr, start, end)
x |
An |
chr |
A |
start |
Unused (for signature compatibility). |
end |
Unused (for signature compatibility). |
A filtered RBedMethyl object.
Summarize methylation by a set of regions.
summarizeByRegion(x, regions)summarizeByRegion(x, regions)
x |
An |
regions |
A |
A DataFrame with coverage, mod_reads, beta, and n_sites.
lines <- c( paste("chr1", 0, 1, "m", 0, "+", 0, 1, 0, 10, 0.5, 5, 5, 0, 0, 0, 0, 0, sep = "\t"), paste("chr1", 10, 11, "m", 0, "+", 10, 11, 0, 20, 0.25, 5, 15, 0, 0, 0, 0, 0, sep = "\t") ) tmp <- tempfile(fileext = ".bed") writeLines(lines, tmp) bm <- readBedMethyl(tmp, mod = "m", fields = c("coverage", "pct", "mod_reads")) regions <- GenomicRanges::GRanges( seqnames = "chr1", ranges = IRanges::IRanges(start = 1, end = 12) ) summarizeByRegion(bm, regions)lines <- c( paste("chr1", 0, 1, "m", 0, "+", 0, 1, 0, 10, 0.5, 5, 5, 0, 0, 0, 0, 0, sep = "\t"), paste("chr1", 10, 11, "m", 0, "+", 10, 11, 0, 20, 0.25, 5, 15, 0, 0, 0, 0, 0, sep = "\t") ) tmp <- tempfile(fileext = ".bed") writeLines(lines, tmp) bm <- readBedMethyl(tmp, mod = "m", fields = c("coverage", "pct", "mod_reads")) regions <- GenomicRanges::GRanges( seqnames = "chr1", ranges = IRanges::IRanges(start = 1, end = 12) ) summarizeByRegion(bm, regions)