1 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.

2 Load the package

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

library(MACSr)

3 Usage

3.1 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.

3.2 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()
#> snapshotDate(): 2021-10-18
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  @ Tue, 26 Oct 2021 17:47:38: 
#> # Command line: 
#> # ARGUMENTS LIST:
#> # name = run_callpeak_narrow0
#> # format = AUTO
#> # ChIP-seq file = ['/home/biocbuild/.cache/R/ExperimentHub/3d73832332a552_4601']
#> # control file = ['/home/biocbuild/.cache/R/ExperimentHub/3d73831c5d77a3_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  @ Tue, 26 Oct 2021 17:47:38: #1 read tag files... 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #1 read treatment tags... 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: Detected format is: BED 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: * Input file is gzipped. 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #1.2 read input tags... 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: Detected format is: BED 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: * Input file is gzipped. 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #1 tag size is determined as 101 bps 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #1 tag size = 101.0 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #1  total tags in treatment: 49622 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #1 user defined the maximum tags... 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #1 filter out redundant tags at the same location and the same strand by allowing at most 1 tag(s) 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #1  tags after filtering in treatment: 48047 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #1  Redundant rate of treatment: 0.03 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #1  total tags in control: 50837 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #1 user defined the maximum tags... 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #1 filter out redundant tags at the same location and the same strand by allowing at most 1 tag(s) 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #1  tags after filtering in control: 50783 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #1  Redundant rate of control: 0.00 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #1 finished! 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #2 Build Peak Model... 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #2 looking for paired plus/minus strand peaks... 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #2 Total number of paired peaks: 469 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #2 Model building with cross-correlation: Done 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #2 finished! 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #2 predicted fragment length is 228 bps 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #2 alternative fragment length(s) may be 228 bps 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #2.2 Generate R script for model : /tmp/RtmpROQspP/run_callpeak_narrow0_model.r 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #3 Call peaks... 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #3 Pre-compute pvalue-qvalue table... 
#> INFO  @ Tue, 26 Oct 2021 17:47:38: #3 Cutoff vs peaks called will be analyzed! 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #3 Analysis of cutoff vs num of peaks or total length has been saved in b'/tmp/RtmpROQspP/run_callpeak_narrow0_cutoff_analysis.txt' 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #3 In the peak calling step, the following will be performed simultaneously: 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #3   Write bedGraph files for treatment pileup (after scaling if necessary)... run_callpeak_narrow0_treat_pileup.bdg 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #3   Write bedGraph files for control lambda (after scaling if necessary)... run_callpeak_narrow0_control_lambda.bdg 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #3   Pileup will be based on sequencing depth in treatment. 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #3 Call peaks for each chromosome... 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #4 Write output xls file... /tmp/RtmpROQspP/run_callpeak_narrow0_peaks.xls 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #4 Write peak in narrowPeak format file... /tmp/RtmpROQspP/run_callpeak_narrow0_peaks.narrowPeak 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #4 Write summits bed file... /tmp/RtmpROQspP/run_callpeak_narrow0_summits.bed 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: Done!
cp2 <- callpeak(CHIP, CTRL, gsize = 5.2e7, store_bdg = TRUE,
                name = "run_callpeak_broad", outdir = tempdir(),
                broad = TRUE)
#> INFO  @ Tue, 26 Oct 2021 17:47:39: 
#> # Command line: 
#> # ARGUMENTS LIST:
#> # name = run_callpeak_broad
#> # format = AUTO
#> # ChIP-seq file = ['/home/biocbuild/.cache/R/ExperimentHub/3d73832332a552_4601']
#> # control file = ['/home/biocbuild/.cache/R/ExperimentHub/3d73831c5d77a3_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  @ Tue, 26 Oct 2021 17:47:39: #1 read tag files... 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #1 read treatment tags... 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: Detected format is: BED 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: * Input file is gzipped. 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #1.2 read input tags... 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: Detected format is: BED 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: * Input file is gzipped. 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #1 tag size is determined as 101 bps 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #1 tag size = 101.0 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #1  total tags in treatment: 49622 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #1 user defined the maximum tags... 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #1 filter out redundant tags at the same location and the same strand by allowing at most 1 tag(s) 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #1  tags after filtering in treatment: 48047 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #1  Redundant rate of treatment: 0.03 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #1  total tags in control: 50837 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #1 user defined the maximum tags... 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #1 filter out redundant tags at the same location and the same strand by allowing at most 1 tag(s) 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #1  tags after filtering in control: 50783 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #1  Redundant rate of control: 0.00 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #1 finished! 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #2 Build Peak Model... 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #2 looking for paired plus/minus strand peaks... 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #2 Total number of paired peaks: 469 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #2 Model building with cross-correlation: Done 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #2 finished! 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #2 predicted fragment length is 228 bps 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #2 alternative fragment length(s) may be 228 bps 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #2.2 Generate R script for model : /tmp/RtmpROQspP/run_callpeak_broad_model.r 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #3 Call peaks... 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #3 Call broad peaks with given level1 -log10qvalue cutoff and level2: 1.301030, 1.000000... 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #3 Pre-compute pvalue-qvalue table... 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #3 In the peak calling step, the following will be performed simultaneously: 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #3   Write bedGraph files for treatment pileup (after scaling if necessary)... run_callpeak_broad_treat_pileup.bdg 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #3   Write bedGraph files for control lambda (after scaling if necessary)... run_callpeak_broad_control_lambda.bdg 
#> INFO  @ Tue, 26 Oct 2021 17:47:39: #3 Call peaks for each chromosome... 
#> INFO  @ Tue, 26 Oct 2021 17:47:40: #4 Write output xls file... /tmp/RtmpROQspP/run_callpeak_broad_peaks.xls 
#> INFO  @ Tue, 26 Oct 2021 17:47:40: #4 Write broad peak in broadPeak format file... /tmp/RtmpROQspP/run_callpeak_broad_peaks.broadPeak 
#> INFO  @ Tue, 26 Oct 2021 17:47:40: #4 Write broad peak in bed12/gappedPeak format file... /tmp/RtmpROQspP/run_callpeak_broad_peaks.gappedPeak 
#> INFO  @ Tue, 26 Oct 2021 17:47:40: Done!

Here are the outputs.

cp1
#> macsList class
#> $outputs:
#>  /tmp/RtmpROQspP/run_callpeak_narrow0_control_lambda.bdg
#>  /tmp/RtmpROQspP/run_callpeak_narrow0_cutoff_analysis.txt
#>  /tmp/RtmpROQspP/run_callpeak_narrow0_model.r
#>  /tmp/RtmpROQspP/run_callpeak_narrow0_peaks.narrowPeak
#>  /tmp/RtmpROQspP/run_callpeak_narrow0_peaks.xls
#>  /tmp/RtmpROQspP/run_callpeak_narrow0_summits.bed
#>  /tmp/RtmpROQspP/run_callpeak_narrow0_treat_pileup.bdg 
#> $arguments: tfile, cfile, gsize, outdir, name, store_bdg, cutoff_analysis 
#> $log:
#>  INFO  @ Tue, 26 Oct 2021 17:47:38: 
#>  # Command line: 
#>  # ARGUMENTS LIST:
#>  # name = run_callpeak_narrow0
#>  # format = AUTO
#> ...
cp2
#> macsList class
#> $outputs:
#>  /tmp/RtmpROQspP/run_callpeak_broad_control_lambda.bdg
#>  /tmp/RtmpROQspP/run_callpeak_broad_model.r
#>  /tmp/RtmpROQspP/run_callpeak_broad_peaks.broadPeak
#>  /tmp/RtmpROQspP/run_callpeak_broad_peaks.gappedPeak
#>  /tmp/RtmpROQspP/run_callpeak_broad_peaks.xls
#>  /tmp/RtmpROQspP/run_callpeak_broad_treat_pileup.bdg 
#> $arguments: tfile, cfile, gsize, outdir, name, store_bdg, broad 
#> $log:
#>  INFO  @ Tue, 26 Oct 2021 17:47:39: 
#>  # Command line: 
#>  # ARGUMENTS LIST:
#>  # name = run_callpeak_broad
#>  # format = AUTO
#> ...

3.3 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/RtmpROQspP/run_callpeak_narrow0_control_lambda.bdg" 
#> [2] "/tmp/RtmpROQspP/run_callpeak_narrow0_cutoff_analysis.txt"
#> [3] "/tmp/RtmpROQspP/run_callpeak_narrow0_model.r"            
#> [4] "/tmp/RtmpROQspP/run_callpeak_narrow0_peaks.narrowPeak"   
#> [5] "/tmp/RtmpROQspP/run_callpeak_narrow0_peaks.xls"          
#> [6] "/tmp/RtmpROQspP/run_callpeak_narrow0_summits.bed"        
#> [7] "/tmp/RtmpROQspP/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(cp1$log)
#> INFO  @ Tue, 26 Oct 2021 17:47:38:  # Command line:  # ARGUMENTS LIST: # name = run_callpeak_narrow0 # format = AUTO # ChIP-seq file = ['/home/biocbuild/.cache/R/ExperimentHub/3d73832332a552_4601'] # control file = ['/home/biocbuild/.cache/R/ExperimentHub/3d73831c5d77a3_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  @ Tue, 26 Oct 2021 17:47:38: #1 read tag files...  INFO  @ Tue, 26 Oct 2021 17:47:38: #1 read treatment tags...  INFO  @ Tue, 26 Oct 2021 17:47:38: Detected format is: BED  INFO  @ Tue, 26 Oct 2021 17:47:38: * Input file is gzipped.  INFO  @ Tue, 26 Oct 2021 17:47:38: #1.2 read input tags...  INFO  @ Tue, 26 Oct 2021 17:47:38: Detected format is: BED  INFO  @ Tue, 26 Oct 2021 17:47:38: * Input file is gzipped.  INFO  @ Tue, 26 Oct 2021 17:47:38: #1 tag size is determined as 101 bps  INFO  @ Tue, 26 Oct 2021 17:47:38: #1 tag size = 101.0  INFO  @ Tue, 26 Oct 2021 17:47:38: #1  total tags in treatment: 49622  INFO  @ Tue, 26 Oct 2021 17:47:38: #1 user defined the maximum tags...  INFO  @ Tue, 26 Oct 2021 17:47:38: #1 filter out redundant tags at the same location and the same strand by allowing at most 1 tag(s)  INFO  @ Tue, 26 Oct 2021 17:47:38: #1  tags after filtering in treatment: 48047  INFO  @ Tue, 26 Oct 2021 17:47:38: #1  Redundant rate of treatment: 0.03  INFO  @ Tue, 26 Oct 2021 17:47:38: #1  total tags in control: 50837  INFO  @ Tue, 26 Oct 2021 17:47:38: #1 user defined the maximum tags...  INFO  @ Tue, 26 Oct 2021 17:47:38: #1 filter out redundant tags at the same location and the same strand by allowing at most 1 tag(s)  INFO  @ Tue, 26 Oct 2021 17:47:38: #1  tags after filtering in control: 50783  INFO  @ Tue, 26 Oct 2021 17:47:38: #1  Redundant rate of control: 0.00  INFO  @ Tue, 26 Oct 2021 17:47:38: #1 finished!  INFO  @ Tue, 26 Oct 2021 17:47:38: #2 Build Peak Model...  INFO  @ Tue, 26 Oct 2021 17:47:38: #2 looking for paired plus/minus strand peaks...  INFO  @ Tue, 26 Oct 2021 17:47:38: #2 Total number of paired peaks: 469  INFO  @ Tue, 26 Oct 2021 17:47:38: #2 Model building with cross-correlation: Done  INFO  @ Tue, 26 Oct 2021 17:47:38: #2 finished!  INFO  @ Tue, 26 Oct 2021 17:47:38: #2 predicted fragment length is 228 bps  INFO  @ Tue, 26 Oct 2021 17:47:38: #2 alternative fragment length(s) may be 228 bps  INFO  @ Tue, 26 Oct 2021 17:47:38: #2.2 Generate R script for model : /tmp/RtmpROQspP/run_callpeak_narrow0_model.r  INFO  @ Tue, 26 Oct 2021 17:47:38: #3 Call peaks...  INFO  @ Tue, 26 Oct 2021 17:47:38: #3 Pre-compute pvalue-qvalue table...  INFO  @ Tue, 26 Oct 2021 17:47:38: #3 Cutoff vs peaks called will be analyzed!  INFO  @ Tue, 26 Oct 2021 17:47:39: #3 Analysis of cutoff vs num of peaks or total length has been saved in b'/tmp/RtmpROQspP/run_callpeak_narrow0_cutoff_analysis.txt'  INFO  @ Tue, 26 Oct 2021 17:47:39: #3 In the peak calling step, the following will be performed simultaneously:  INFO  @ Tue, 26 Oct 2021 17:47:39: #3   Write bedGraph files for treatment pileup (after scaling if necessary)... run_callpeak_narrow0_treat_pileup.bdg  INFO  @ Tue, 26 Oct 2021 17:47:39: #3   Write bedGraph files for control lambda (after scaling if necessary)... run_callpeak_narrow0_control_lambda.bdg  INFO  @ Tue, 26 Oct 2021 17:47:39: #3   Pileup will be based on sequencing depth in treatment.  INFO  @ Tue, 26 Oct 2021 17:47:39: #3 Call peaks for each chromosome...  INFO  @ Tue, 26 Oct 2021 17:47:39: #4 Write output xls file... /tmp/RtmpROQspP/run_callpeak_narrow0_peaks.xls  INFO  @ Tue, 26 Oct 2021 17:47:39: #4 Write peak in narrowPeak format file... /tmp/RtmpROQspP/run_callpeak_narrow0_peaks.narrowPeak  INFO  @ Tue, 26 Oct 2021 17:47:39: #4 Write summits bed file... /tmp/RtmpROQspP/run_callpeak_narrow0_summits.bed  INFO  @ Tue, 26 Oct 2021 17:47:39: Done!

4 Resources

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

5 SessionInfo

sessionInfo()
#> R version 4.1.1 (2021-08-10)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 20.04.3 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.14-bioc/R/lib/libRblas.so
#> LAPACK: /home/biocbuild/bbs-3.14-bioc/R/lib/libRlapack.so
#> 
#> locale:
#>  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
#>  [3] LC_TIME=en_GB              LC_COLLATE=C              
#>  [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       
#> 
#> attached base packages:
#> [1] stats     graphics  grDevices utils     datasets  methods   base     
#> 
#> other attached packages:
#> [1] MACSdata_1.1.0   MACSr_1.2.0      BiocStyle_2.22.0
#> 
#> loaded via a namespace (and not attached):
#>  [1] Biobase_2.54.0                httr_1.4.2                   
#>  [3] sass_0.4.0                    bit64_4.0.5                  
#>  [5] jsonlite_1.7.2                AnnotationHub_3.2.0          
#>  [7] here_1.0.1                    bslib_0.3.1                  
#>  [9] shiny_1.7.1                   assertthat_0.2.1             
#> [11] interactiveDisplayBase_1.32.0 BiocManager_1.30.16          
#> [13] stats4_4.1.1                  BiocFileCache_2.2.0          
#> [15] blob_1.2.2                    GenomeInfoDbData_1.2.7       
#> [17] yaml_2.2.1                    BiocVersion_3.14.0           
#> [19] pillar_1.6.4                  RSQLite_2.2.8                
#> [21] lattice_0.20-45               glue_1.4.2                   
#> [23] reticulate_1.22               digest_0.6.28                
#> [25] XVector_0.34.0                promises_1.2.0.1             
#> [27] htmltools_0.5.2               httpuv_1.6.3                 
#> [29] Matrix_1.3-4                  pkgconfig_2.0.3              
#> [31] dir.expiry_1.2.0              zlibbioc_1.40.0              
#> [33] bookdown_0.24                 purrr_0.3.4                  
#> [35] xtable_1.8-4                  later_1.3.0                  
#> [37] tibble_3.1.5                  KEGGREST_1.34.0              
#> [39] generics_0.1.1                IRanges_2.28.0               
#> [41] ellipsis_0.3.2                withr_2.4.2                  
#> [43] cachem_1.0.6                  BiocGenerics_0.40.0          
#> [45] magrittr_2.0.1                crayon_1.4.1                 
#> [47] mime_0.12                     memoise_2.0.0                
#> [49] evaluate_0.14                 fansi_0.5.0                  
#> [51] tools_4.1.1                   lifecycle_1.0.1              
#> [53] basilisk.utils_1.6.0          stringr_1.4.0                
#> [55] S4Vectors_0.32.0              AnnotationDbi_1.56.0         
#> [57] Biostrings_2.62.0             compiler_4.1.1               
#> [59] jquerylib_0.1.4               GenomeInfoDb_1.30.0          
#> [61] rlang_0.4.12                  grid_4.1.1                   
#> [63] RCurl_1.98-1.5                rappdirs_0.3.3               
#> [65] bitops_1.0-7                  rmarkdown_2.11               
#> [67] basilisk_1.6.0                ExperimentHub_2.2.0          
#> [69] DBI_1.1.1                     curl_4.3.2                   
#> [71] R6_2.5.1                      knitr_1.36                   
#> [73] dplyr_1.0.7                   fastmap_1.1.0                
#> [75] bit_4.0.4                     utf8_1.2.2                   
#> [77] rprojroot_2.0.2               filelock_1.0.2               
#> [79] stringi_1.7.5                 parallel_4.1.1               
#> [81] Rcpp_1.0.7                    vctrs_0.3.8                  
#> [83] png_0.1-7                     dbplyr_2.1.1                 
#> [85] tidyselect_1.1.1              xfun_0.27