Introduction

With the improvement of sequencing techniques, chromatin immunoprecipitation followed by high throughput sequencing (ChIP-Seq) is getting popular to study genome-wide protein-DNA interactions. To address the lack of powerful ChIP-Seq analysis method, we presented the Model-based Analysis of ChIP-Seq (MACS), for identifying transcript factor binding sites. MACS captures the influence of genome complexity to evaluate the significance of enriched ChIP regions and MACS improves the spatial resolution of binding sites through combining the information of both sequencing tag position and orientation. MACS can be easily used for ChIP-Seq data alone, or with a control sample with the increase of specificity. Moreover, as a general peak-caller, MACS can also be applied to any “DNA enrichment assays” if the question to be asked is simply: where we can find significant reads coverage than the random background.

This package is a wrapper of the MACS toolkit based on basilisk.

Load the package

The package is built on basilisk. The dependent python library macs3 will be installed automatically inside its conda environment.

library(MACSr)

Usage

MACS3 functions

There are 13 functions imported from MACS3. Details of each function can be checked from its manual.

Functions Description
callpeak Main MACS3 Function to call peaks from alignment results.
bdgpeakcall Call peaks from bedGraph output.
bdgbroadcall Call broad peaks from bedGraph output.
bdgcmp Comparing two signal tracks in bedGraph format.
bdgopt Operate the score column of bedGraph file.
cmbreps Combine BEDGraphs of scores from replicates.
bdgdiff Differential peak detection based on paired four bedGraph files.
filterdup Remove duplicate reads, then save in BED/BEDPE format.
predictd Predict d or fragment size from alignment results.
pileup Pileup aligned reads (single-end) or fragments (paired-end)
randsample Randomly choose a number/percentage of total reads.
refinepeak Take raw reads alignment, refine peak summits.
callvar Call variants in given peak regions from the alignment BAM files.
hmmratac Dedicated peak calling based on Hidden Markov Model for ATAC-seq data.

Function callpeak

We have uploaded multipe test datasets from MACS to a data package MACSdata in the ExperimentHub. For example, Here we download a pair of single-end bed files to run the callpeak function.

eh <- ExperimentHub::ExperimentHub()
eh <- AnnotationHub::query(eh, "MACSdata")
CHIP <- eh[["EH4558"]]
#> see ?MACSdata and browseVignettes('MACSdata') for documentation
#> loading from cache
CTRL <- eh[["EH4563"]]
#> see ?MACSdata and browseVignettes('MACSdata') for documentation
#> loading from cache

Here is an example to call narrow and broad peaks on the SE bed files.

cp1 <- callpeak(CHIP, CTRL, gsize = 5.2e7, store_bdg = TRUE,
                name = "run_callpeak_narrow0", outdir = tempdir(),
                cutoff_analysis = TRUE)
#> INFO  @ 17 Nov 2023 14:29:22: [635 MB] 
#> # Command line: 
#> # ARGUMENTS LIST:
#> # name = run_callpeak_narrow0
#> # format = AUTO
#> # ChIP-seq file = ['/home/qhu/.cache/R/ExperimentHub/1071e4528bfa6_4601']
#> # control file = ['/home/qhu/.cache/R/ExperimentHub/1071e71fe08e5_4606']
#> # effective genome size = 5.20e+07
#> # band width = 300
#> # model fold = [5.0, 50.0]
#> # qvalue cutoff = 5.00e-02
#> # The maximum gap between significant sites is assigned as the read length/tag size.
#> # The minimum length of peaks is assigned as the predicted fragment length "d".
#> # Larger dataset will be scaled towards smaller dataset.
#> # Range for calculating regional lambda is: 1000 bps and 10000 bps
#> # Broad region calling is off
#> # Additional cutoff on fold-enrichment is: 0.10
#> # Paired-End mode is off
#>  
#> INFO  @ 17 Nov 2023 14:29:22: [635 MB] #1 read tag files... 
#> INFO  @ 17 Nov 2023 14:29:22: [635 MB] #1 read treatment tags... 
#> INFO  @ 17 Nov 2023 14:29:22: [638 MB] Detected format is: BED 
#> INFO  @ 17 Nov 2023 14:29:22: [638 MB] * Input file is gzipped. 
#> INFO  @ 17 Nov 2023 14:29:22: [642 MB] #1.2 read input tags... 
#> INFO  @ 17 Nov 2023 14:29:22: [642 MB] Detected format is: BED 
#> INFO  @ 17 Nov 2023 14:29:22: [642 MB] * Input file is gzipped. 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1 tag size is determined as 101 bps 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1 tag size = 101.0 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1  total tags in treatment: 49622 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1 user defined the maximum tags... 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1 filter out redundant tags at the same location and the same strand by allowing at most 1 tag(s) 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1  tags after filtering in treatment: 48047 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1  Redundant rate of treatment: 0.03 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1  total tags in control: 50837 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1 user defined the maximum tags... 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1 filter out redundant tags at the same location and the same strand by allowing at most 1 tag(s) 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1  tags after filtering in control: 50783 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1  Redundant rate of control: 0.00 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1 finished! 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #2 Build Peak Model... 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #2 looking for paired plus/minus strand peaks... 
#> INFO  @ 17 Nov 2023 14:29:23: [643 MB] #2 Total number of paired peaks: 469 
#> INFO  @ 17 Nov 2023 14:29:23: [643 MB] #2 Model building with cross-correlation: Done 
#> INFO  @ 17 Nov 2023 14:29:23: [643 MB] #2 finished! 
#> INFO  @ 17 Nov 2023 14:29:23: [643 MB] #2 predicted fragment length is 228 bps 
#> INFO  @ 17 Nov 2023 14:29:23: [643 MB] #2 alternative fragment length(s) may be 228 bps 
#> INFO  @ 17 Nov 2023 14:29:23: [643 MB] #2.2 Generate R script for model : /tmp/RtmpIUsM2o/run_callpeak_narrow0_model.r 
#> INFO  @ 17 Nov 2023 14:29:23: [643 MB] #3 Call peaks... 
#> INFO  @ 17 Nov 2023 14:29:23: [643 MB] #3 Pre-compute pvalue-qvalue table... 
#> INFO  @ 17 Nov 2023 14:29:23: [643 MB] #3 Cutoff vs peaks called will be analyzed! 
#> INFO  @ 17 Nov 2023 14:29:23: [661 MB] #3 Analysis of cutoff vs num of peaks or total length has been saved in b'/tmp/RtmpIUsM2o/run_callpeak_narrow0_cutoff_analysis.txt' 
#> INFO  @ 17 Nov 2023 14:29:23: [661 MB] #3 In the peak calling step, the following will be performed simultaneously: 
#> INFO  @ 17 Nov 2023 14:29:23: [661 MB] #3   Write bedGraph files for treatment pileup (after scaling if necessary)... run_callpeak_narrow0_treat_pileup.bdg 
#> INFO  @ 17 Nov 2023 14:29:23: [661 MB] #3   Write bedGraph files for control lambda (after scaling if necessary)... run_callpeak_narrow0_control_lambda.bdg 
#> INFO  @ 17 Nov 2023 14:29:23: [661 MB] #3   Pileup will be based on sequencing depth in treatment. 
#> INFO  @ 17 Nov 2023 14:29:23: [661 MB] #3 Call peaks for each chromosome... 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #4 Write output xls file... /tmp/RtmpIUsM2o/run_callpeak_narrow0_peaks.xls 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #4 Write peak in narrowPeak format file... /tmp/RtmpIUsM2o/run_callpeak_narrow0_peaks.narrowPeak 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #4 Write summits bed file... /tmp/RtmpIUsM2o/run_callpeak_narrow0_summits.bed 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] Done!
cp2 <- callpeak(CHIP, CTRL, gsize = 5.2e7, store_bdg = TRUE,
                name = "run_callpeak_broad", outdir = tempdir(),
                broad = TRUE)
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] 
#> # Command line: 
#> # ARGUMENTS LIST:
#> # name = run_callpeak_broad
#> # format = AUTO
#> # ChIP-seq file = ['/home/qhu/.cache/R/ExperimentHub/1071e4528bfa6_4601']
#> # control file = ['/home/qhu/.cache/R/ExperimentHub/1071e71fe08e5_4606']
#> # effective genome size = 5.20e+07
#> # band width = 300
#> # model fold = [5.0, 50.0]
#> # qvalue cutoff for narrow/strong regions = 5.00e-02
#> # qvalue cutoff for broad/weak regions = 1.00e-01
#> # The maximum gap between significant sites is assigned as the read length/tag size.
#> # The minimum length of peaks is assigned as the predicted fragment length "d".
#> # Larger dataset will be scaled towards smaller dataset.
#> # Range for calculating regional lambda is: 1000 bps and 10000 bps
#> # Broad region calling is on
#> # Additional cutoff on fold-enrichment is: 0.10
#> # Paired-End mode is off
#>  
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #1 read tag files... 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #1 read treatment tags... 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] Detected format is: BED 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] * Input file is gzipped. 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #1.2 read input tags... 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] Detected format is: BED 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] * Input file is gzipped. 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #1 tag size is determined as 101 bps 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #1 tag size = 101.0 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #1  total tags in treatment: 49622 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #1 user defined the maximum tags... 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #1 filter out redundant tags at the same location and the same strand by allowing at most 1 tag(s) 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #1  tags after filtering in treatment: 48047 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #1  Redundant rate of treatment: 0.03 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #1  total tags in control: 50837 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #1 user defined the maximum tags... 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #1 filter out redundant tags at the same location and the same strand by allowing at most 1 tag(s) 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #1  tags after filtering in control: 50783 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #1  Redundant rate of control: 0.00 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #1 finished! 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #2 Build Peak Model... 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #2 looking for paired plus/minus strand peaks... 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #2 Total number of paired peaks: 469 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #2 Model building with cross-correlation: Done 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #2 finished! 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #2 predicted fragment length is 228 bps 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #2 alternative fragment length(s) may be 228 bps 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #2.2 Generate R script for model : /tmp/RtmpIUsM2o/run_callpeak_broad_model.r 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #3 Call peaks... 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #3 Call broad peaks with given level1 -log10qvalue cutoff and level2: 1.301030, 1.000000... 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #3 Pre-compute pvalue-qvalue table... 
#> INFO  @ 17 Nov 2023 14:29:24: [663 MB] #3 In the peak calling step, the following will be performed simultaneously: 
#> INFO  @ 17 Nov 2023 14:29:24: [663 MB] #3   Write bedGraph files for treatment pileup (after scaling if necessary)... run_callpeak_broad_treat_pileup.bdg 
#> INFO  @ 17 Nov 2023 14:29:24: [663 MB] #3   Write bedGraph files for control lambda (after scaling if necessary)... run_callpeak_broad_control_lambda.bdg 
#> INFO  @ 17 Nov 2023 14:29:24: [663 MB] #3 Call peaks for each chromosome... 
#> INFO  @ 17 Nov 2023 14:29:24: [664 MB] #4 Write output xls file... /tmp/RtmpIUsM2o/run_callpeak_broad_peaks.xls 
#> INFO  @ 17 Nov 2023 14:29:24: [664 MB] #4 Write broad peak in broadPeak format file... /tmp/RtmpIUsM2o/run_callpeak_broad_peaks.broadPeak 
#> INFO  @ 17 Nov 2023 14:29:24: [664 MB] #4 Write broad peak in bed12/gappedPeak format file... /tmp/RtmpIUsM2o/run_callpeak_broad_peaks.gappedPeak 
#> INFO  @ 17 Nov 2023 14:29:24: [664 MB] Done!

Here are the outputs.

cp1
#> macsList class
#> $outputs:
#>  /tmp/RtmpIUsM2o/run_callpeak_narrow0_control_lambda.bdg
#>  /tmp/RtmpIUsM2o/run_callpeak_narrow0_cutoff_analysis.txt
#>  /tmp/RtmpIUsM2o/run_callpeak_narrow0_model.r
#>  /tmp/RtmpIUsM2o/run_callpeak_narrow0_peaks.narrowPeak
#>  /tmp/RtmpIUsM2o/run_callpeak_narrow0_peaks.xls
#>  /tmp/RtmpIUsM2o/run_callpeak_narrow0_summits.bed
#>  /tmp/RtmpIUsM2o/run_callpeak_narrow0_treat_pileup.bdg 
#> $arguments: tfile, cfile, gsize, outdir, name, store_bdg, cutoff_analysis 
#> $log:
#>  INFO  @ 17 Nov 2023 14:29:22: [635 MB] 
#>  # Command line: 
#>  # ARGUMENTS LIST:
#>  # name = run_callpeak_narrow0
#>  # format = AUTO
#> ...
cp2
#> macsList class
#> $outputs:
#>  /tmp/RtmpIUsM2o/run_callpeak_broad_control_lambda.bdg
#>  /tmp/RtmpIUsM2o/run_callpeak_broad_model.r
#>  /tmp/RtmpIUsM2o/run_callpeak_broad_peaks.broadPeak
#>  /tmp/RtmpIUsM2o/run_callpeak_broad_peaks.gappedPeak
#>  /tmp/RtmpIUsM2o/run_callpeak_broad_peaks.xls
#>  /tmp/RtmpIUsM2o/run_callpeak_broad_treat_pileup.bdg 
#> $arguments: tfile, cfile, gsize, outdir, name, store_bdg, broad 
#> $log:
#>  INFO  @ 17 Nov 2023 14:29:23: [662 MB] 
#>  # Command line: 
#>  # ARGUMENTS LIST:
#>  # name = run_callpeak_broad
#>  # format = AUTO
#> ...

The macsList class

The macsList is designed to contain everything of an execution, including function, inputs, outputs and logs, for the purpose of reproducibility.

For example, we can the function and input arguments.

cp1$arguments
#> [[1]]
#> callpeak
#> 
#> $tfile
#> CHIP
#> 
#> $cfile
#> CTRL
#> 
#> $gsize
#> [1] 5.2e+07
#> 
#> $outdir
#> tempdir()
#> 
#> $name
#> [1] "run_callpeak_narrow0"
#> 
#> $store_bdg
#> [1] TRUE
#> 
#> $cutoff_analysis
#> [1] TRUE

The files of all the outputs are collected.

cp1$outputs
#> [1] "/tmp/RtmpIUsM2o/run_callpeak_narrow0_control_lambda.bdg" 
#> [2] "/tmp/RtmpIUsM2o/run_callpeak_narrow0_cutoff_analysis.txt"
#> [3] "/tmp/RtmpIUsM2o/run_callpeak_narrow0_model.r"            
#> [4] "/tmp/RtmpIUsM2o/run_callpeak_narrow0_peaks.narrowPeak"   
#> [5] "/tmp/RtmpIUsM2o/run_callpeak_narrow0_peaks.xls"          
#> [6] "/tmp/RtmpIUsM2o/run_callpeak_narrow0_summits.bed"        
#> [7] "/tmp/RtmpIUsM2o/run_callpeak_narrow0_treat_pileup.bdg"

The log is especially important for MACS to check. Detailed information was given in the log when running.

cat(paste(cp1$log, collapse="\n"))
#> INFO  @ 17 Nov 2023 14:29:22: [635 MB] 
#> # Command line: 
#> # ARGUMENTS LIST:
#> # name = run_callpeak_narrow0
#> # format = AUTO
#> # ChIP-seq file = ['/home/qhu/.cache/R/ExperimentHub/1071e4528bfa6_4601']
#> # control file = ['/home/qhu/.cache/R/ExperimentHub/1071e71fe08e5_4606']
#> # effective genome size = 5.20e+07
#> # band width = 300
#> # model fold = [5.0, 50.0]
#> # qvalue cutoff = 5.00e-02
#> # The maximum gap between significant sites is assigned as the read length/tag size.
#> # The minimum length of peaks is assigned as the predicted fragment length "d".
#> # Larger dataset will be scaled towards smaller dataset.
#> # Range for calculating regional lambda is: 1000 bps and 10000 bps
#> # Broad region calling is off
#> # Additional cutoff on fold-enrichment is: 0.10
#> # Paired-End mode is off
#>  
#> INFO  @ 17 Nov 2023 14:29:22: [635 MB] #1 read tag files... 
#> INFO  @ 17 Nov 2023 14:29:22: [635 MB] #1 read treatment tags... 
#> INFO  @ 17 Nov 2023 14:29:22: [638 MB] Detected format is: BED 
#> INFO  @ 17 Nov 2023 14:29:22: [638 MB] * Input file is gzipped. 
#> INFO  @ 17 Nov 2023 14:29:22: [642 MB] #1.2 read input tags... 
#> INFO  @ 17 Nov 2023 14:29:22: [642 MB] Detected format is: BED 
#> INFO  @ 17 Nov 2023 14:29:22: [642 MB] * Input file is gzipped. 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1 tag size is determined as 101 bps 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1 tag size = 101.0 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1  total tags in treatment: 49622 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1 user defined the maximum tags... 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1 filter out redundant tags at the same location and the same strand by allowing at most 1 tag(s) 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1  tags after filtering in treatment: 48047 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1  Redundant rate of treatment: 0.03 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1  total tags in control: 50837 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1 user defined the maximum tags... 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1 filter out redundant tags at the same location and the same strand by allowing at most 1 tag(s) 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1  tags after filtering in control: 50783 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1  Redundant rate of control: 0.00 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #1 finished! 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #2 Build Peak Model... 
#> INFO  @ 17 Nov 2023 14:29:22: [643 MB] #2 looking for paired plus/minus strand peaks... 
#> INFO  @ 17 Nov 2023 14:29:23: [643 MB] #2 Total number of paired peaks: 469 
#> INFO  @ 17 Nov 2023 14:29:23: [643 MB] #2 Model building with cross-correlation: Done 
#> INFO  @ 17 Nov 2023 14:29:23: [643 MB] #2 finished! 
#> INFO  @ 17 Nov 2023 14:29:23: [643 MB] #2 predicted fragment length is 228 bps 
#> INFO  @ 17 Nov 2023 14:29:23: [643 MB] #2 alternative fragment length(s) may be 228 bps 
#> INFO  @ 17 Nov 2023 14:29:23: [643 MB] #2.2 Generate R script for model : /tmp/RtmpIUsM2o/run_callpeak_narrow0_model.r 
#> INFO  @ 17 Nov 2023 14:29:23: [643 MB] #3 Call peaks... 
#> INFO  @ 17 Nov 2023 14:29:23: [643 MB] #3 Pre-compute pvalue-qvalue table... 
#> INFO  @ 17 Nov 2023 14:29:23: [643 MB] #3 Cutoff vs peaks called will be analyzed! 
#> INFO  @ 17 Nov 2023 14:29:23: [661 MB] #3 Analysis of cutoff vs num of peaks or total length has been saved in b'/tmp/RtmpIUsM2o/run_callpeak_narrow0_cutoff_analysis.txt' 
#> INFO  @ 17 Nov 2023 14:29:23: [661 MB] #3 In the peak calling step, the following will be performed simultaneously: 
#> INFO  @ 17 Nov 2023 14:29:23: [661 MB] #3   Write bedGraph files for treatment pileup (after scaling if necessary)... run_callpeak_narrow0_treat_pileup.bdg 
#> INFO  @ 17 Nov 2023 14:29:23: [661 MB] #3   Write bedGraph files for control lambda (after scaling if necessary)... run_callpeak_narrow0_control_lambda.bdg 
#> INFO  @ 17 Nov 2023 14:29:23: [661 MB] #3   Pileup will be based on sequencing depth in treatment. 
#> INFO  @ 17 Nov 2023 14:29:23: [661 MB] #3 Call peaks for each chromosome... 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #4 Write output xls file... /tmp/RtmpIUsM2o/run_callpeak_narrow0_peaks.xls 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #4 Write peak in narrowPeak format file... /tmp/RtmpIUsM2o/run_callpeak_narrow0_peaks.narrowPeak 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] #4 Write summits bed file... /tmp/RtmpIUsM2o/run_callpeak_narrow0_summits.bed 
#> INFO  @ 17 Nov 2023 14:29:23: [662 MB] Done!

Resources

More details about MACS3 can be found: https://macs3-project.github.io/MACS/.

SessionInfo

sessionInfo()
#> R version 4.3.1 (2023-06-16)
#> Platform: x86_64-conda-linux-gnu (64-bit)
#> Running under: Ubuntu 22.04.3 LTS
#> 
#> Matrix products: default
#> BLAS/LAPACK: /home/qhu/miniconda3/envs/rcwl/lib/libopenblasp-r0.3.23.so;  LAPACK version 3.11.0
#> 
#> locale:
#>  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
#>  [3] LC_TIME=en_US.UTF-8        LC_COLLATE=en_US.UTF-8    
#>  [5] LC_MONETARY=en_US.UTF-8    LC_MESSAGES=en_US.UTF-8   
#>  [7] LC_PAPER=en_US.UTF-8       LC_NAME=C                 
#>  [9] LC_ADDRESS=C               LC_TELEPHONE=C            
#> [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C       
#> 
#> time zone: America/New_York
#> tzcode source: system (glibc)
#> 
#> attached base packages:
#> [1] stats     graphics  grDevices utils     datasets  methods   base     
#> 
#> other attached packages:
#> [1] MACSdata_1.8.0   MACSr_1.11.2     BiocStyle_2.28.1
#> 
#> loaded via a namespace (and not attached):
#>  [1] tidyselect_1.2.0              dplyr_1.1.2                  
#>  [3] blob_1.2.4                    bitops_1.0-7                 
#>  [5] filelock_1.0.2                Biostrings_2.68.1            
#>  [7] fastmap_1.1.1                 RCurl_1.98-1.13              
#>  [9] BiocFileCache_2.8.0           promises_1.2.0.1             
#> [11] digest_0.6.33                 mime_0.12                    
#> [13] lifecycle_1.0.3               ellipsis_0.3.2               
#> [15] KEGGREST_1.40.1               interactiveDisplayBase_1.40.0
#> [17] RSQLite_2.3.1                 magrittr_2.0.3               
#> [19] compiler_4.3.1                rlang_1.1.1                  
#> [21] sass_0.4.7                    tools_4.3.1                  
#> [23] utf8_1.2.3                    yaml_2.3.7                   
#> [25] knitr_1.43                    bit_4.0.5                    
#> [27] curl_5.0.1                    here_1.0.1                   
#> [29] reticulate_1.30               withr_2.5.0                  
#> [31] purrr_1.0.1                   BiocGenerics_0.46.0          
#> [33] desc_1.4.2                    grid_4.3.1                   
#> [35] stats4_4.3.1                  fansi_1.0.4                  
#> [37] ExperimentHub_2.10.0          xtable_1.8-4                 
#> [39] cli_3.6.1                     rmarkdown_2.25               
#> [41] crayon_1.5.2                  ragg_1.2.5                   
#> [43] generics_0.1.3                httr_1.4.6                   
#> [45] DBI_1.1.3                     cachem_1.0.8                 
#> [47] stringr_1.5.0                 zlibbioc_1.46.0              
#> [49] parallel_4.3.1                AnnotationDbi_1.62.2         
#> [51] BiocManager_1.30.21.1         XVector_0.40.0               
#> [53] basilisk_1.12.1               vctrs_0.6.3                  
#> [55] Matrix_1.6-0                  jsonlite_1.8.7               
#> [57] dir.expiry_1.8.0              bookdown_0.36                
#> [59] IRanges_2.34.1                S4Vectors_0.38.1             
#> [61] bit64_4.0.5                   systemfonts_1.0.4            
#> [63] jquerylib_0.1.4               glue_1.6.2                   
#> [65] pkgdown_2.0.7                 stringi_1.7.12               
#> [67] BiocVersion_3.17.1            later_1.3.1                  
#> [69] GenomeInfoDb_1.36.4           tibble_3.2.1                 
#> [71] pillar_1.9.0                  basilisk.utils_1.12.1        
#> [73] rappdirs_0.3.3                htmltools_0.5.5              
#> [75] GenomeInfoDbData_1.2.10       R6_2.5.1                     
#> [77] dbplyr_2.3.3                  textshaping_0.3.6            
#> [79] rprojroot_2.0.3               evaluate_0.21                
#> [81] shiny_1.7.4.1                 lattice_0.21-8               
#> [83] Biobase_2.60.0                AnnotationHub_3.10.0         
#> [85] png_0.1-8                     memoise_2.0.1                
#> [87] httpuv_1.6.11                 bslib_0.5.0                  
#> [89] Rcpp_1.0.11                   xfun_0.39                    
#> [91] fs_1.6.3                      pkgconfig_2.0.3