Title: | MeRIP-Seq data Analysis for Genomic Power Investigation and Evaluation |
---|---|
Description: | This package aims to perform power analysis for the MeRIP-seq study. It calculates FDR, FDC, power, and precision under various study design parameters, including but not limited to sample size, sequencing depth, and testing method. It can also output results into .xlsx files or produce corresponding figures of choice. |
Authors: | Daoyu Duan [aut, cre], Zhenxing Guo [aut] |
Maintainer: | Daoyu Duan <[email protected]> |
License: | MIT + file LICENSE |
Version: | 1.7.0 |
Built: | 2024-11-18 03:41:57 UTC |
Source: | https://github.com/bioc/magpie |
This function plots all power measurements of a certain sequencing depth in a 2x2 panel. Power measurements to plot include "FDR", "FDC", "Power", and "Precision".
plotAll(Power.list, depth_factor = 1)
plotAll(Power.list, depth_factor = 1)
Power.list |
A list produced by |
depth_factor |
A numerical value indicating which sequencing depth to plot. For example, 2 means doubling the original sequencing depth. Default is 1. |
It plots all power measurements of a certain sequencing depth in a 2x2 panel. Power measurements to plot include "FDR", "FDC", "Power", and "Precision".
library(magpie) ### Main function power.test <- quickPower(dataset = "GSE46705", test_method = "TRESS") ### plot all in a panel under sequencing depth 1x plotAll(power.test, depth_factor = 1)
library(magpie) ### Main function power.test <- quickPower(dataset = "GSE46705", test_method = "TRESS") ### plot all in a panel under sequencing depth 1x plotAll(power.test, depth_factor = 1)
This function plots all power measurements of the original sequencing depth by strata in a 2x2 panel. Power measurements to plot include "FDR", "FDC", "Power", and "Precision".
plotAll_Strata(Power.list)
plotAll_Strata(Power.list)
Power.list |
A list produced by |
It plots all power measurements of the original sequencing depth by strata in a 2x2 panel. Power measurements to plot include "FDR", "FDC", "Power", and "Precision".
library(magpie) ### Main function power.test <- quickPower(dataset = "GSE55575", test_method = "TRESS") ### plot all strata results in a panel plotAll_Strata(power.test)
library(magpie) ### Main function power.test <- quickPower(dataset = "GSE55575", test_method = "TRESS") ### plot all strata results in a panel plotAll_Strata(power.test)
This function plots a certain power measurement of a certain sequencing depth. Power measurements to plot include "FDR", "FDC", "Power", and "Precision".
plotRes(Power.list, depth_factor = 1, value_option = "FDR")
plotRes(Power.list, depth_factor = 1, value_option = "FDR")
Power.list |
A list produced by |
depth_factor |
A numerical value indicating which sequencing depth to plot. For example, 2 means doubling the original sequencing depth. Default is 1. |
value_option |
A character indicating which measurement to plot. Options include "FDR", "FDC", "Power", and "Precision". |
It plots a certain power measurement of a certain sequencing depth. Power measurements to plot include "FDR", "FDC", "Power", and "Precision".
library(magpie) ### Main function power.test <- quickPower(dataset = "GSE46705", test_method = "TRESS") ### plot FDR under sequencing depth 1x plotRes(power.test, depth_factor = 1, value_option = "FDR")
library(magpie) ### Main function power.test <- quickPower(dataset = "GSE46705", test_method = "TRESS") ### plot FDR under sequencing depth 1x plotRes(power.test, depth_factor = 1, value_option = "FDR")
This function plots a certain power measurement of the original sequencing depth by strata. Power measurements to plot include "FDR", "FDC", "Power", and "Precision".
plotStrata(Power.list, value_option = "FDR")
plotStrata(Power.list, value_option = "FDR")
Power.list |
A list produced by |
value_option |
A character indicating which measurement to plot. Options include "FDR", "FDC", "Power", and "Precision". |
It plots a certain power measurement of the original sequencing depth by strata. Power measurements to plot include "FDR", "FDC", "Power", and "Precision".
library(magpie) ### Main function power.test <- quickPower(dataset = "GSE46705", test_method = "TRESS") ### plot a FDR strata result plotStrata(power.test, value_option = "FDR")
library(magpie) ### Main function power.test <- quickPower(dataset = "GSE46705", test_method = "TRESS") ### plot a FDR strata result plotStrata(power.test, value_option = "FDR")
Human HeLa cell line: Two replicates of wild type (WT) and two replicates of knockdown (KD) of complex METTL3.
data(GSE46705_TRESS_res)
data(GSE46705_TRESS_res)
A list object.
data(GSE46705_TRESS_res)
data(GSE46705_TRESS_res)
This function conducts simulations with various user-defined study design parameters, including but not limited to sample size, sequencing depth, and testing method. Users will need to provide either partial or whole-genome MeRIP-seq data for parameter estimation purposes.
powerEval( Input.file, IP.file, BamDir, annoDir, variable, bam_factor, nsim = 10, N.reps = c(2, 3, 5), depth_factor = c(1, 2, 5), thres = c(0.01, 0.05, 0.1, 0.2), Test_method = "TRESS" )
powerEval( Input.file, IP.file, BamDir, annoDir, variable, bam_factor, nsim = 10, N.reps = c(2, 3, 5), depth_factor = c(1, 2, 5), thres = c(0.01, 0.05, 0.1, 0.2), Test_method = "TRESS" )
Input.file |
A vector containing the name of BAM files of input samples. |
IP.file |
A vector containing the name of BAM files of IP samples. |
BamDir |
A character stating the directory path of all .BAM files. |
annoDir |
A character stating the directory path of the ".sqlite" file for annotation. |
variable |
A vector indicating the experimental conditions of all samples. |
bam_factor |
A numrical value indicating the ratio of provided data to the whole genome data. Default is 0.05. |
nsim |
An integer indicating the number of iterations to simulate under each scenario. Default is 10. |
N.reps |
A vector of integers indicating the numbers of replicates to simulate, in both groups. Default is c(2,3,5). |
depth_factor |
A vector of numerical values indicating how much sequencing depth of the provided data will be increased in simulations. For example, 2 means doubling the original sequencing depth. Default is c(1,2,5). |
thres |
A vector of numerical values indicating the p-value thresholds used in power calculation. Default is c(0.01, 0.05, 0.1, 0.2). |
Test_method |
A character indicating which DMR calling method to use. Options are "TRESS" and "exomePeak2". Default is "TRESS". |
A list of calculated power measurements that will be used as the input of functions writeToxlsx
, writeToxlsx_strata
, plotAll
, plotRes
, plotAll_Strata
, and plotStrata
.
Measurements include:
FDR |
The ratio of number of false positives to the number of positive discoveries. |
FDC |
The ratio of number of false positives to the number of true positives. |
Power |
Statistical power. |
Precision |
The ratio of number of true positives to the number of positive discoveries. |
## Not run: library(magpieData) library(magpie) ### Get the example data BAM_path <- getBAMpath() ### Call PowerEval() power.test <- powerEval( Input.file = c("Ctrl1.chr15.input.bam", "Ctrl2.chr15.input.bam", "Case1.chr15.input.bam", "Case2.chr15.input.bam"), IP.file = c("Ctrl1.chr15.ip.bam", "Ctrl2.chr15.ip.bam", "Case1.chr15.ip.bam", "Case2.chr15.ip.bam"), BamDir = BAM_path, annoDir = paste0(BAM_path, "/hg18_chr15.sqlite"), variable = rep(c("Ctrl", "Trt"), each = 2), bam_factor = 0.03, nsim = 10, N.reps = c(2, 3, 5, 7), depth_factor = c(1, 2), thres = c(0.01, 0.05, 0.1), Test_method = "TRESS" ) ## End(Not run)
## Not run: library(magpieData) library(magpie) ### Get the example data BAM_path <- getBAMpath() ### Call PowerEval() power.test <- powerEval( Input.file = c("Ctrl1.chr15.input.bam", "Ctrl2.chr15.input.bam", "Case1.chr15.input.bam", "Case2.chr15.input.bam"), IP.file = c("Ctrl1.chr15.ip.bam", "Ctrl2.chr15.ip.bam", "Case1.chr15.ip.bam", "Case2.chr15.ip.bam"), BamDir = BAM_path, annoDir = paste0(BAM_path, "/hg18_chr15.sqlite"), variable = rep(c("Ctrl", "Trt"), each = 2), bam_factor = 0.03, nsim = 10, N.reps = c(2, 3, 5, 7), depth_factor = c(1, 2), thres = c(0.01, 0.05, 0.1), Test_method = "TRESS" ) ## End(Not run)
This function quickly outputs pre-calculated power evaluation results from four GEO MeRIP-seq datasets: (GSE46705, GSE55575, GSE115105, and GSE94613). The obtained results can be used to generate Excel files and various figures.
quickPower(dataset = "GSE46705", test_method = "TRESS")
quickPower(dataset = "GSE46705", test_method = "TRESS")
dataset |
A character specifying the selected dataset. Default is 'GSE46705'. Options are 'GSE46705', 'GSE55575', 'GSE115105', and 'GSE94613'. |
test_method |
A character indicating which DMR calling method to use. Options are "TRESS" and "exomePeak2". Default is "TRESS". |
GSE46705: Human HeLa cell line: Two replicates of wild type (WT) and two replicates of knockdown (KD) of complex METTL3.
GSE55575: Mouse embryonic fibroblasts: Two replicates of wild type (WT) and four replicates of knockdown (KD) of WTAP.
GSE115105: Two sample types from WT and YTHDF1 KO mice. Each type has two replicates.
GSE94613: Human leukemia cell line: Four replicates of wild type (WT) and eight replicates of knockdown (KD) of complex METTL3.
A list of calculated power measurements that will be used as the input of functions writeToxlsx
, writeToxlsx_strata
, plotAll
, plotRes
, plotAll_Strata
, and plotStrata
.
Measurements include:
FDR |
The ratio of number of false positives to the number of positive discoveries. |
FDC |
The ratio of number of false positives to the number of true positives. |
Power |
Statistical power. |
Precision |
The ratio of number of true positives to the number of positive discoveries. |
library(magpie) power.test <- quickPower(dataset = "GSE46705")
library(magpie) power.test <- quickPower(dataset = "GSE46705")
This function writes power evaulation results to a .xlsx file.
writeToxlsx(pl, file)
writeToxlsx(pl, file)
pl |
A list produced by |
file |
A character indicating the name of the output .xlsx file. |
It outputs a .xlsx file including FDR, FDC, power, and precision under various sample sizes, sequencing depths, and adjusted p-value thresholds.
library(magpie) ### Main function power.test <- quickPower(dataset = "GSE46705", test_method = "TRESS") ### write out .xlsx writeToxlsx(power.test, file = "test_TRESS.xlsx")
library(magpie) ### Main function power.test <- quickPower(dataset = "GSE46705", test_method = "TRESS") ### write out .xlsx writeToxlsx(power.test, file = "test_TRESS.xlsx")
This function writes power evaulation results of four strata to a .xlsx file. Only results from the original sequencing depth are saved. Here, strata are determined by mean input control levels of simulated data.
writeToxlsx_strata(pl, file)
writeToxlsx_strata(pl, file)
pl |
A list produced by |
file |
A character indicating the name of the output .xlsx file. |
It outputs a .xlsx file including FDR, FDC, power, and precision under the original sequencing depth and various sample sizes and input stratas.
library(magpie) ### Main function power.test <- quickPower(dataset = "GSE46705", test_method = "TRESS") ### write out .xlsx writeToxlsx_strata(power.test, file = "test_strata_TRESS.xlsx")
library(magpie) ### Main function power.test <- quickPower(dataset = "GSE46705", test_method = "TRESS") ### write out .xlsx writeToxlsx_strata(power.test, file = "test_strata_TRESS.xlsx")