Analyzes the enrichment of ranges within R-loop forming sequences (RLFS). See details.
analyzeRLFS(
object,
mask = NULL,
quiet = FALSE,
useMask = TRUE,
noZ = FALSE,
ntimes = 100,
stepsize = 50,
...
)
An RLRanges object.
GRanges object containing masked genomic ranges. Not needed unless masked genome unavailable (see genomeMasks). Custom masks can be generated using regioneR::getMask.
If TRUE, messages are suppressed. Default: FALSE.
If FALSE, masked genome is not used. This is not recommended unless a mask is unavailable as it can lead to spurious results. Default: TRUE.
If TRUE, Z-score distribution is not calculated. Default: FALSE.
Number of permutations to perform (default: 100).
The step size for calculating the Z score distribution. Default: 50. See also regioneR::localZScore.
Arguments passed to regioneR::permTest.
An RLRanges object with RLFS analysis results
accessible via RLSeq::rlresult(object, "rlfsRes")
. Contains the
following structure:
perTestResults
An object of the class permTestResultsList
from regioneR
with
the results of permutation testing. See also
regioneR::permTest for full description.
Z-scores
An object of the class localZScoreResultsList
from regioneR
.
Contains the results of local Z-score analysis +/-5kb around each RLFS.
See also regioneR::localZScore.
R-loop forming sequences are regions of the genome with sequences that are favorable for R-loop formation. They are computationally predicted with the QmRLFS-finder software program and serve as a data-independent test of whether a sample has mapped R-loops robustly or not.
Permutation testing is implemented via regioneR::permTest such that, for each permutation, R-loop peaks were randomized using regioneR::circularRandomizeRegions and then the number of overlaps with RLFS are counted. 100 permutations are used by default to build an empirical distribution for peak/RLFS overlap. Then the true number of overlaps from non-randomized peaks and RLFS are compared to the null distribution to calculate Z-score and significance of enrichment. Finally, a Z-score distribution was calculated (using regioneR::localZScore) 5kb upstream and downstream of the average RLFS midpoint.
These results are subsequently used in the binary classification of the sample as "POS" (maps R-loops) or "NEG" (does not map R-loops). See also predictCondition.
# Example dataset
rlr <- readRDS(system.file("extdata", "rlrsmall.rds", package = "RLSeq"))
# Perform RLFS analysis (remove ntimes=2 and noZ=TRUE for a typical analysis)
rlr <- analyzeRLFS(rlr, ntimes = 2, noZ = TRUE)
#> - Evaluating Inputs...
#> - Running permTest...
#> [1] "Note: The minimum p-value with only 2 permutations is 0.333333333333333. You should consider increasing the number of permutations."
#> Warning: All permuted values are equal to 1. Z-score is infinite.
#> - Extracting pileup...
#> - Done.