Description Usage Arguments Details Value See Also Examples
scanMeripBAM check and organize the BAM files in MeRIP-seq data before peak calling using exomePeakCalling.
The library types of the RNA-seq and the filters such as SAM FLAG score are specified in this function.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20  | scanMeripBAM(
  bam_ip = NULL,
  bam_input = NULL,
  bam_treated_ip = NULL,
  bam_treated_input = NULL,
  paired_end = FALSE,
  library_type = c("unstranded", "1st_strand", "2nd_strand"),
  index_bam = TRUE,
  bam_files = NULL,
  design_ip = NULL,
  design_treatment = NULL,
  mapq = 30L,
  isSecondaryAlignment = FALSE,
  isNotPassingQualityControls = FALSE,
  isDuplicate = FALSE,
  isPaired = NA,
  isProperPair = NA,
  hasUnmappedMate = NA,
  ...
)
 | 
bam_ip | 
 a   | 
bam_input | 
 a   | 
bam_treated_ip | 
 a   | 
bam_treated_input | 
 a  If the bam files do not contain treatment group, user should only fill the arguments of   | 
paired_end | 
 a   | 
library_type | 
 a  
  | 
index_bam | 
 a  The BAM index files will be named by adding ".bai" after the names of the corresponding BAM files.  | 
bam_files | 
 optional, a   | 
design_ip | 
 optional, a   | 
design_treatment | 
 optional, a   | 
mapq | 
 a non-negative integer specifying the minimum reads mapping quality. BAM records with mapping qualities less than   | 
isPaired, isProperPair, hasUnmappedMate, isSecondaryAlignment, isNotPassingQualityControls, isDuplicate, ... | 
 arguments specifying the filters on SAM FLAG scores, inherited from   | 
scanMeripBAM takes the input of the BAM file directories for the MeRIP-seq datasets.
It first checks the completeness of the BAM files and the BAM indexes. Then, the design information of IP/input and treated/control are returned as a MeripBamFileList object.
If the BAM file indexes are missing, the BAM files will be automatically indexed with the package Rsamtools.
a MeripBamFileList object.
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  | ### Define BAM File Directories
f1 = system.file("extdata", "IP1.bam", package="exomePeak2")
f2 = system.file("extdata", "IP2.bam", package="exomePeak2")
f3 = system.file("extdata", "IP3.bam", package="exomePeak2")
f4 = system.file("extdata", "IP4.bam", package="exomePeak2")
IP_BAM = c(f1,f2,f3,f4)
f1 = system.file("extdata", "Input1.bam", package="exomePeak2")
f2 = system.file("extdata", "Input2.bam", package="exomePeak2")
f3 = system.file("extdata", "Input3.bam", package="exomePeak2")
INPUT_BAM = c(f1,f2,f3)
f1 = system.file("extdata", "treated_IP1.bam", package="exomePeak2")
TREATED_IP_BAM = c(f1)
f1 = system.file("extdata", "treated_Input1.bam", package="exomePeak2")
TREATED_INPUT_BAM = c(f1)
### For MeRIP-Seq Experiment Without the Treatment Group
MeRIP_Seq_Alignment <- scanMeripBAM(
  bam_ip = IP_BAM,
  bam_input = INPUT_BAM,
  paired_end = FALSE
)
### For MeRIP-Seq Experiment With the Treatment Group
MeRIP_Seq_Alignment <- scanMeripBAM(
bam_ip = IP_BAM,
bam_input = INPUT_BAM,
bam_treated_ip = TREATED_IP_BAM,
bam_treated_input = TREATED_INPUT_BAM,
paired_end = FALSE
)
 | 
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.