Description Usage Arguments Details Value Note Author(s) See Also Examples
View source: R/logScoreWithControl.R
Compute the RIPScore using both RIP and control posteriors for each bins, merge and summarize the scores for the merged bins, and finally compute the p-value and adjusted p-value for the summary RIPScore.
1 | logScoreWithControl(nbhGRRIP, nbhGRCTL, padjMethod = "BH", getControlStats = TRUE)
|
nbhGRRIP |
|
nbhGRCTL |
An optional arugment as a |
padjMethod |
Method used to adjust multiple testing performed in |
getControlStats |
Binary indicator to whether return statistics including computed for the control library alone. If TRUE, then all control specific score columns will be reported with a prefix "CTL". |
The RIPScore is compupted in computeLogOdd
as the log odd ratio of the posterior for the RIP state (z_{i} = 2) over the posterior for the background state (z_{i} = 1) in RIP library subtracted by the log odd ratio computed from the control library. The adjacent bins with hidden states predicted by nbh_vit
as the enriched state (corresponding to the NB with larger mean) are merged. The RIPSscores are averaged over the merged bins. To assess the statistical significance of the RIPScore for each region, we assume that the RIPScore follows a Gaussian (Normal) distribution with mean and standard deviation estimated using the RIPScores over all of the bins. The rationale is based on the assumption that most of the RIPScores correspond to the background state and together contribute to a stable estimate of the test statistics (TS) and p-value computed using the R built-in function pnorm
. The p-value is adjusted by p.adjust
with BH method by default. The same procedure is applied optionally to the control library.
GRanges
of merged bins with values slot saved for RIPScore (lodOdd), p-value (pval), adjusted p-value (pvalAdj) for RIP and optionally for control.
Internal function used by seekRIP
.
Yue Li
logScoreWithoutControl, seekRIP, computeLogOdd, scoreMergedBins
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 | if(interactive()) { # check the example in seekRIP
# Retrieve system files
extdata.dir <- system.file("extdata", package="RIPSeeker")
bamFiles <- list.files(extdata.dir, ".bam$", recursive=TRUE, full.names=TRUE)
bamFiles <- grep("PRC2", bamFiles, value=TRUE)
# Parameters setting
binSize <- 1e5 # use a large fixed bin size for demo only
multicore <- FALSE # use multicore
strandType <- "-" # set strand type to minus strand
################ run main function for HMM inference on all chromosomes ################
mainSeekOutputRIP <-
mainSeek(bamFiles=grep(pattern="SRR039214",
bamFiles, value=TRUE, invert=TRUE),
binSize=binSize, strandType=strandType,
reverseComplement=TRUE, genomeBuild="mm9",
uniqueHit = TRUE, assignMultihits = TRUE,
rerunWithDisambiguatedMultihits = FALSE,
multicore=multicore, silentMain=FALSE, verbose=TRUE)
mainSeekOutputCTL <- mainSeek(bamFiles=grep(pattern="SRR039214",
bamFiles, value=TRUE, invert=FALSE),
binSize=binSize, strandType=strandType,
reverseComplement=TRUE, genomeBuild="mm9",
uniqueHit = TRUE, assignMultihits = TRUE,
rerunWithDisambiguatedMultihits = FALSE,
multicore=multicore, silentMain=FALSE, verbose=TRUE)
################ Compute log score and test for significance WITH control ################
ripGR.wicontrol <- logScoreWithControl(mainSeekOutputRIP$nbhGRList$chrX, mainSeekOutputCTL$nbhGRList$chrX)
ripGR.wicontrol
}
|
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.