Quantify windows by summing across positions in each of Nmod and Nvalid
Source: R/scanGenome.R
sumNmodNvalid.RdCalculate sums of per-position Nmod and Nvalid counts for
positions that overlap windows in gr, and calculate FracMod
as the ratio between the two.
Arguments
- se
A
RangedSummarizedExperimentobject, as generated byreadModBam(..., level = "summary")(i.e. with the assays"Nmod"and"Nvalid").- gr
A
GRangesobject defining the windows to quantify.- includeEmpty
Logical scalar. If
TRUE, include also windows without overlapping modified bases in the output (with zero counts in both the Nmod and Nvalid assays).
Value
A RangedSummarizedExperiment
object with one row per window defined by gr and the window
quantifications contained in the assays "Nmod", "Nvalid"
and "FracMod". This object is suitable for differential analysis
by getDifferentiallyModifiedWindows.
Examples
library(GenomicRanges)
library(SummarizedExperiment)
modbamfiles <- system.file("extdata",
c("6mA_1_10reads.bam", "6mA_2_10reads.bam"),
package = "footprintR")
se <- readModBam(bamfiles = modbamfiles, regions = "chr1:6940000-6940500",
level = "summary", modbase = "a", trim = TRUE,
BPPARAM = BiocParallel::SerialParam())
s <- seq(6940000, 6940400, by = 100)
windowgr <- GRanges(seqnames = "chr1",
ranges = IRanges(start = s, width = 100, names = letters[1:5]))
seSum <- sumNmodNvalid(se = se, gr = windowgr)
seSum
#> class: RangedSummarizedExperiment
#> dim: 5 2
#> metadata(1): readLevelData
#> assays(3): Nmod Nvalid FracMod
#> rownames(5): a b c d e
#> rowData names(0):
#> colnames(2): s1 s2
#> colData names(2): sample modbase
assay(seSum, "Nvalid")
#> s1 s2
#> a 99 55
#> b 91 63
#> c 51 52
#> d 82 68
#> e 76 41