Quantify windows by calculating a Fourier-transform-based phasing score
Source:R/scanGenome.R
phasingScoreFourier.Rd
Calculate a phasing score by decomposing FracMod
in a window into
frequency components and extract the numCoef
-th coefficient.
Arguments
- se
A
RangedSummarizedExperiment
object, as generated byreadModBam(..., level = "summary")
(i.e. with the assays"Nvalid"
and"FracMod"
).- gr
A
GRanges
object defining the windows to quantify. Here, the width of the ranges ingr
have to be an integer multiple of the period of interest (typicallywidth(gr) = (numCoef - 1) * periodOfInterest
).- numCoef
An integer scalar giving the Fourier coefficient to extract. The coefficient corresponds to a period of
width(gr) / (numCoef - 1)
.
Value
A RangedSummarizedExperiment
object with one row per window defined by gr
and the window
phasing scores for periodic signals with a frequency defined by
"numCoef"
.
This object is suitable for transformation by
getRangesWithAssayValues
.
Examples
library(GenomicRanges)
library(SummarizedExperiment)
modbamfiles <- system.file("extdata",
c("6mA_1_10reads.bam", "6mA_2_10reads.bam"),
package = "footprintR")
rng <- GRanges("chr1", IRanges(6940000, 6943600))
windowStep <- 180
windowSize <- 4 * windowStep
s <- seq(start(rng), end(rng) - windowSize + 1, by = windowStep)
windowgr <- GRanges(seqnames = seqnames(rng),
ranges = IRanges(start = s, width = windowSize))
se <- readModBam(bamfiles = modbamfiles, regions = rng, level = "summary",
modbase = "a", trim = TRUE,
BPPARAM = BiocParallel::SerialParam())
seFourier <- phasingScoreFourier(se = se, gr = windowgr, numCoef = 5)
assay(seFourier, "phasingScoreAbs") # amplitude of requested coefficient
#> s1 s2
#> [1,] 22.104124 11.7576205
#> [2,] 27.132726 3.8888579
#> [3,] 17.332390 0.4179011
#> [4,] 9.322762 0.0000000
#> [5,] 8.247887 0.0000000
#> [6,] 4.598294 0.0000000
#> [7,] 0.984381 0.0000000
#> [8,] 0.000000 0.0000000
#> [9,] 0.000000 0.0000000
#> [10,] 0.000000 0.0000000
#> [11,] 0.000000 0.0000000
#> [12,] 0.000000 0.0000000
#> [13,] 0.000000 0.0000000
#> [14,] 0.000000 0.0000000
#> [15,] 0.000000 0.0000000
#> [16,] 0.000000 0.0000000
#> [17,] 0.000000 0.0000000
assay(seFourier, "phasingScoreRel") # relative amplitude
#> s1 s2
#> [1,] 0.032935952 0.008596142
#> [2,] 0.027635236 0.005619616
#> [3,] 0.015046715 0.007232441
#> [4,] 0.009746856 NaN
#> [5,] 0.010512054 NaN
#> [6,] 0.008950576 NaN
#> [7,] 0.006229053 NaN
#> [8,] NaN NaN
#> [9,] NaN NaN
#> [10,] NaN NaN
#> [11,] NaN NaN
#> [12,] NaN NaN
#> [13,] NaN NaN
#> [14,] NaN NaN
#> [15,] NaN NaN
#> [16,] NaN NaN
#> [17,] NaN NaN