The esATAC package provides a dataflow graphs organized pipeline for quantifying and annotating ATAC-seq and DNase-seq Reads in R, which integrate the functionality of several R packages (such as Rsamtools, ChIPpeakAnno and so on) and external softwares (e.g. AdapterRemoval[1], bowtie2[2], through the Rowtie2 package and Fseq[3]). Users can call the preset pipeline functions or rebuild the workflow with elements process of esATAC for sequence data (FASTQ, SAM, BAM, BED file) complete/partial processing and analyzing easily and flexibly in a single R script. That will be convenient to migrate, share and reproduce all details such as parameters settings, intermediate result and so on. Besides, a quality control report file in HTML, which is able to be viewed in web browser, will be created in preset pipelines.
esATAC can be easily installed on various operator system platforms (Windows, Linux, Mac OS). All functions in package consume up to 16G memory. Most function only consume less than 8G. So the package is available for not only servers but also most of PC. esATAC supports analysis of both single end reads and paired-end ATAC-seq and DNase-seq sequencing reads data. We have successfully apply the pipeline on raw datasets (FASTQ files) from GEO. The related sequencing platform are Illumina. All valid standard format intermediate result files (FASTQ, SAM, BAM, BED file) generated by other program (such as BAM BED files from ENCODE) are also tested by rebuilt sub-pipeline.
If you do not know where to start with ATAC-seq or DNase-seq data, you can print flowchart like this:
library(esATAC)
printMap()
Following the flowchart, related functions could be found in manual. For example,if you want to query “SamToBed” in the flowchart, you can add “atac” prefix to query mannual like this:
?atacSamToBed
or use the lowercase of the initial letter
?samToBed
The workflow start with “UnzipAndMerge” function atacUnzipAndMerge. It unzips and merges the replicates into one FASTQ file(two for paired end file). Names of reads will be renamed as numbers: 1,2,3,… by calling “Renamer” function atacRenamer. The file will be smaller for further analysis. Adapter of reads may be found and removed by “RemoveAdapter” function atacRemoveAdapter. Then reads are ready for mapping to reference genome. “Bowtie2Mapping” mapping function atacBowtie2Mapping can do this job. “SamToBam”, “Rsortbam”,“BamToBed”,“SamToBed” and “BedUtils” provide general processing methods for SAM file including converting format into BAM or BED file, sorting according to chromosome/start site/end site, reads conditional filtering, reads shifting and so on. The ready-use reads in BED file may call peak by “PeakCallingFseq” function atacPeakCallingFseq. The pipeline also provide quality control elements (e.g. “FastQC”, “LibComplexQC”, “TSSQC”) and some general genome function analysis elements (e.g. “RMotifScan”,“RPeakAnno”). For more detail, you can see the manual or the examples in following sections.
This package is developed and maintained by members of Xiaowo Wang Lab
MOE Key Laboratory of Bioinformatics and Bioinformatics Division,
TNLIST / Department of Automation, Tsinghua University
email:{wei-z14,w-zhang16}(at)mails.tsinghua.edu.cn
To install the latest version of esATAC, you will need to be using the latest version of R. esATAC is part of Bioconductor project, so you can install esATAC and its dependencies like this:
source("http://www.bioconductor.org/biocLite.R")
biocLite("esATAC")
Just like other R package, you need to load esATAC like this each time before using the package.
library(esATAC)
If you need to use fseq, we recommend to set max memory size for java (8G, 8000M in the example). Or rJava will use the default parameter for fseq.
options(java.parameters = "-Xmx8000m")
The BSgenome package, TxDb known gene package and OrgDb annotation package for some functions are required. We recommend to install and load the specific species related packages before using the packages.
library(magrittr)
library(BSgenome.Hsapiens.UCSC.hg19)
library(TxDb.Hsapiens.UCSC.hg19.knownGene)
library(org.Hs.eg.db)
library(R.utils)
These configurations are also optional. “tmpdir” is the path to save all of the temporary data and the default result storage path. If it is not configured, current work directory will be set as “tmpdir”. “threads” is the maximum threads allowed to be created for data processing. The default value is 1. More thread will consume more memery in some processes.
# we use temp directiory "td" here
# Change it to your directiory because the intermediate file may be huge
td<-tempdir()
options(atacConf=setConfigure("tmpdir",td))
options(atacConf=setConfigure("threads",8))
We strongly recommend to install reference data first before using the package although it is optional. “refdir” is the folder that will save all of the reference data. “genome” is the genome name like hg19, hg38, mm10, mm9 and so on. The program will detect the elements that have not been installed and install them. Some resources need to be downloaded from internet. So don’t forget to connect internet during installation. Or the installation will be failed. If all of the reference data was installed, these two lines still need to be called for configuring the reference data path and genome.
#uncomment and modify to run:
#options(atacConf=setConfigure("refdir","path/to/refdatafolder"))
#options(atacConf=setConfigure("genome","hg19"))
NOTE: The installation will consume several hours for data download and building bowtie2 index depending on computer performance and network bandwidth.
NOTE: The installation is network based. Please keep your network connection. But you don need to worry about disconnect. The program will continue to check finished part and only build unfinish part.
Genome hg19 and mm10 are supported currently. More genomes will be supported in the future.
WARNNING: If the reference data is not configured, the related reference argument of functions has to be set manually during using.
Most of the test datasets in this package (esATAC/extdata/) are generated from GEO: SRR891271 (from GSE47753)[7]. The data is ATAC-seq paired end sequencing for GM12878 cell line. We random sampling 20000 mapped fragments from chr20 and rebuild raw paired-end FASTQ files(file names with chr20 prefix). We also subsample the reads in SAM file and peak calling BED result files. Besides, the files in “uzmg”, “adrm” and “bt2” are the test files from AdapterRemoval and Bowtie2. For detail, you can read the subsequent sections.
esATAC provides an easy-to-use entry, you only need to provide your ATAC-seq sequencing files (fastq format), and assign the spaces and genome assembly, it will do everything for you.
# set reference installation path
options(atacConf=setConfigure("refdir","your/reference/path"))
# set reference genome like hg19, mm10 etc.
options(atacConf=setConfigure("genome","hg19"))
# fastq files path
caseFastqs <- c(
system.file(package="esATAC", "extdata", "chr20_1.1.fq.gz"),
system.file(package="esATAC", "extdata", "chr20_2.1.fq.gz")
)
ctrlFastqs <- c(
system.file(package="esATAC", "extdata", "chr20_1.2.fq.bz2"),
system.file(package="esATAC", "extdata", "chr20_2.2.fq.bz2")
)
# call the esATAC pipeline
atacPipe2(case=list(fastqInput1 = case[1],fastqInput2 = case[2]),
control=list(fastqInput1 = ctrl[1],fastqInput2 = ctrl[2]))
# set reference installation path
options(atacConf=setConfigure("refdir","your/reference/path"))
# set reference genome like hg19, mm10 etc.
options(atacConf=setConfigure("genome","hg19"))
# fastq files path
mate1Fastq <- system.file(package="esATAC", "extdata", "chr20_1.1.fq.gz")
mate2Fastq <- system.file(package="esATAC", "extdata", "chr20_1.2.fq.bz2")
# call the esATAC pipeline
conclusion <- atacPipe(fastqInput1 = mate1Fastq, fastqInput2 = mate2Fastq)
Example report esATAC will download the genome sequence and annotation files, build bowtie2 index, mapping the reads, do the quality control analysis, find peak regions, perform GO analysis and motif enrichment analysis, etc. automatically. Finally, you will get an report file in html format includes to the analysis results.
Download genome sequences and build bowtie index may take some time. If you already have all or some of these files, you can let esATAC skip these steps by renaming them following the format (genome+suffix) and put them in reference installation path.
Example: hg19 reference files
NOTE: genome.1.bt2, genome.2.bt2, genome.3.bt2, genome.4.bt2, genome.rev.1.bt2, genome.rev.2.bt2 are bowtie2 index files; genome.fa is genome sequences. You can put them into reference installation path mannually. genome.DHS.bed, genome.knownGene.sqlite, genome.blacklist.bed are small files, which will be automatically generated or download.
By default, esATAC will perform footprint analysis for all the TF motif PWM matrix in JASPAR database. This step may take a few hours to 2-days for human genome analysis depends on your hardware. If you only want to analyze a specific motif or your own PWMs, you can simple do it like this
motiflist<-list("the pwms")
atacPipe2(case=list(fastqInput1 = case[1],fastqInput2 = case[2]),
control=list(fastqInput1 = ctrl[1],fastqInput2 = ctrl[2]),motifPWM = motiflist)
conclusion <- atacPipe(fastqInput1 = mate1Fastq, fastqInput2 = mate2Fastq,motifPWM = motiflist)
To run the whole pipeline for a typical ATACseq data set of human sample may take ~2 days on a personal computer with single thread. If your thread or R has been stopped during the process. You can simplely resume the analysis by retypin your command line. The program will automatically check the steps that have been finished and continue the anlalysis.
All sub-processes are available for recombine new whole pipeline or sub-pipeline easily and flexibly. They are also able to be called individually. We just show some functions and their combinations from the package. For detail, the users can read the manual.
Users can use %>% to build a pipeline to obtain merged, renamed and adapter removed clean reads fastq file(s) that is ready for mapping.
# Identify adapters
prefix<-system.file(package="esATAC", "extdata", "uzmg")
(reads_1 <-file.path(prefix,"m1",dir(file.path(prefix,"m1"))))
(reads_2 <-file.path(prefix,"m2",dir(file.path(prefix,"m2"))))
reads_merged_1 <- file.path(td,"reads1.fastq")
reads_merged_2 <- file.path(td,"reads2.fastq")
atacproc <-
atacUnzipAndMerge(fastqInput1 = reads_1,fastqInput2 = reads_2) %>%
atacRenamer %>% atacRemoveAdapter
If you want to modify the parameters of AdapterRemoval, you have to refer to Rbowtie2 package:
library(Rbowtie2)
adapterremoval_usage()
If the reference has not been configured, the bowtie2 index should be built first. Then bowtie2 mapping functions could used to map reads to reference genome.
## Building a bowtie2 index
library("Rbowtie2")
refs <- dir(system.file(package="esATAC", "extdata", "bt2","refs"),
full=TRUE)
bowtie2_build(references=refs, bt2Index=file.path(td, "lambda_virus"),
"--threads 4 --quiet",overwrite=TRUE)
## Alignments
reads_1 <- system.file(package="esATAC", "extdata", "bt2", "reads",
"reads_1.fastq")
reads_2 <- system.file(package="esATAC", "extdata", "bt2", "reads",
"reads_2.fastq")
if(file.exists(file.path(td, "lambda_virus.1.bt2"))){
(bowtie2Mapping(bt2Idx = file.path(td, "lambda_virus"),
samOutput = file.path(td, "result.sam"),
fastqInput1=reads_1,fastqInput2=reads_2,threads=3))
head(readLines(file.path(td, "result.sam")))
}
If you want to modify the parameters of bowtie2, you have to refer to Rbowtie2 package:
library(Rbowtie2)
bowtie2_usage()
The mapping results are stored in a SAM file. SamToBed functions can covert it into BED file. During converting, the operation like sorting, shifting, filtering chromosome and so on can also be set to do in the meantime.
sambzfile <- system.file(package="esATAC", "extdata", "Example.sam.bz2")
samfile <- file.path(td,"Example.sam")
bunzip2(sambzfile,destname=samfile,overwrite=TRUE,remove=FALSE)
samToBed(samInput = samfile)
Filter the nucleosome free reads(<100bp) for peak calling.
bedbzfile <- system.file(package="esATAC", "extdata", "chr20.50000.bed.bz2")
bedfile <- file.path(td,"chr20.50000.bed")
bunzip2(bedbzfile,destname=bedfile,overwrite=TRUE,remove=FALSE)
bedUtils(bedInput = bedfile,maxFragLen = 100, chrFilterList = NULL) %>%
atacPeakCalling
ATAC-seq peak locate at open chromatin regions. Annotating these peak could find whether they locate at functional regions(such as promoter and enhancer).
Function “atacPeakAnno” and “peakanno” use function “annotatePeak” in package “ChIPseeker” to annotate ATAC-seq peak. for more information about package “ChIPseeker”, please click here[4].
Function “atacPeakAnno” and “peakanno” accept a bed file path as an input, users can change the parameters like “tssRegion”, “TxDb” according to their require. Now, bioconductor offers many species’ annotation database, click here to search more.
The following example is to exhibit how to annotate a UCSC bed file.
## extract example peak file from package "esATAC"
p1bz <- system.file("extdata", "Example_peak1.bed.bz2", package="esATAC")
peak1_path <- as.vector(bunzip2(filename = p1bz,
destname = file.path(getwd(), "Example_peak1.bed"),
ext="bz2", FUN=bzfile, overwrite=TRUE, remove = FALSE))
## run peakanno to annotate peaks
AnnoInfo <- peakanno(peakInput = peak1_path, TxDb = TxDb.Hsapiens.UCSC.hg19.knownGene, annoDb = "org.Hs.eg.db")
The output contains a pie chart in pdf format like below. It reports the percentage of peaks located in different functional regions.
The function also generate a file(with suffix .df) contains annotation for all peaks. It is converted from dataframe in R, and users could open it with text editor or excel. Below is a part of the output.
chromatin | start | end | annotation | geneStart | geneEnd | geneId | distanceToTSS | SYMBOL |
---|---|---|---|---|---|---|---|---|
chr1 | 416606 | 416895 | Distal Intergenic | 367659 | 368597 | 729759 | 48947 | OR4F29 |
chr1 | 2313275 | 2313587 | Intron (uc001ajb.1/79906, intron 6 of 13) | 2252696 | 2322993 | 79906 | 9406 | MORN1 |
chr1 | 2516858 | 2518618 | Promoter | 2517899 | 2522908 | 127281 | 0 | FAM213B |
chr1 | 2685755 | 2686581 | Intron (uc021oey.1/100287898, intron 4 of 6) | 2572807 | 2706230 | 100287898 | 19649 | TTC34 |
chr1 | 3418588 | 3418822 | Intron (uc001akk.3/1953, intron 14 of 29) | 3404506 | 3448012 | 1953 | 29190 | MEGF6 |
GO analysis is performing enrichment analysis on gene sets. It establishes the relationship between gene sets and functions, and report the most significant function to users.
Function “atacGOAnalysis” and “goanalysis” use function “enrichGO” in package “clusterProfiler” to do GO analysis. for more information about package “clusterProfiler”, please click here[5].
The function need gene Id set as input. User could choose different GO terms(molecular function, biological process and cellular component) according to different input of parameter “ont”.
The following example is to exhibit how to do GO analysis on a gene set.
## extract gene ID
library(clusterProfiler)
data(geneList)
geneId <- names(geneList)[1:100]
## do GO analysis
goAna <- goanalysis(gene = geneId, OrgDb = "org.Hs.eg.db", keytype = "ENTREZID", ont = "MF")
The output file(suffix .df) contains the GO term sorted by p-value, below is a part of the output.
ID | Description | GeneRatio | pvalue | qvalue |
---|---|---|---|---|
GO:0008017 | microtubule binding | 13/93 | 0.0e+00 | 0.00e+00 |
GO:0015631 | tubulin binding | 13/93 | 0.0e+00 | 7.00e-07 |
GO:0003777 | microtubule motor activity | 7/93 | 2.0e-07 | 1.54e-05 |
GO:0050786 | RAGE receptor binding | 4/93 | 3.0e-07 | 1.54e-05 |
GO:0045236 | CXCR chemokine receptor binding | 4/93 | 1.5e-06 | 6.63e-05 |
This function search motif occurrence in the given regions.
Function “atacMotifScan” and “motifscan” use function “matchPWM” in package “Biostrings”, for more parameters and usage, click here[6].
Multi-motif is supported, and the output file is named by your input PWM list. for Multi-motif, we offer parallel computing method for accelerating. Users could specify the parameter “n.cores” to accelerate the program.
The input motif PWM matrix is stored in a list like below.
pwm <- readRDS(system.file("extdata", "motifPWM.rds", package="esATAC"))
pwm
## $CTCF
## 1 2 3 4 5 6
## A 0.02960540 0.03743389 0.04368501 0.02430735 0.001023094 0.05537836
## C 0.04410295 0.03573708 0.02271820 0.05624947 0.057703951 0.00675723
## G 0.02797499 0.04833945 0.04931320 0.01254142 -0.069438277 0.02610485
## T 0.04957887 0.03879208 0.03478530 0.01903153 -0.015525860 0.03014738
## 7 8 9 10 11 12
## A 0.020271900 0.03209952 0.057005153 -0.00460166 0.045748514 0.023909107
## C 0.051259116 0.04889070 0.004838716 -0.06943828 -0.010700941 0.005881734
## G 0.045758368 0.02246973 0.017613999 0.05773064 0.052120515 0.050726179
## T 0.004838716 0.04540858 0.010683413 -0.01070094 0.002433989 0.046034249
## 13 14 15 16 17 18
## A 0.005881734 0.024346457 0.03179891 0.04710491 0.02895828 0.03323194
## C -0.069438277 0.001023094 0.05525243 0.00684137 0.05022974 0.04538490
## G 0.057569617 0.055907445 -0.00460166 0.05082564 0.04481264 0.02756395
## T 0.001023094 0.027215451 0.02651920 0.01005852 0.01942021 0.04786956
## 19
## A 0.04804997
## C 0.03846515
## G 0.04309319
## T 0.02501100
##
## $ATF3
## 1 2 3 4 5 6
## A 0.04367299 0.03800285 0.07142085 0.009661673 -0.02889619 0.07216227
## C 0.05222528 0.03098095 0.03725448 -0.028896186 -0.02889619 -0.02889619
## G 0.04093099 0.07155594 0.04294621 -0.028896186 0.07216227 -0.02889619
## T 0.06989254 0.03510088 0.01473202 0.072148910 -0.02889619 -0.02889619
## 7 8 9 10 11 12
## A 0.01301289 0.02882151 -0.02889619 0.02959694 0.07211734 0.02431670
## C 0.07214318 -0.02889619 -0.02889619 0.07195987 -0.02889619 0.04391574
## G -0.02889619 0.07206064 0.01651982 0.01907820 -0.02889619 0.04041291
## T -0.02889619 -0.02889619 0.07213457 0.02292425 0.02109342 0.07121981
## 13 14
## A 0.03937329 0.06758346
## C 0.07143835 0.04184396
## G 0.03039048 0.05912350
## T 0.03823980 0.04864436
Using “motifscan” function to search motif in given genome regions, UCSC bed file is recommended.
sample.path <- system.file("extdata", "chr20_sample_peak.bed.bz2", package="esATAC")
sample.path <- as.vector(bunzip2(filename = sample.path,
destname = file.path(getwd(), "chr20_sample_peak.bed"),
ext="bz2", FUN=bzfile, overwrite=TRUE, remove = FALSE))
motif.data <- motifscan(peak = sample.path, genome = BSgenome.Hsapiens.UCSC.hg19,
motifPWM = pwm, prefix = "test")
This function reports the exact motif position in the given genome like below(motif: CTCF).
chromatin | start | end | strand | score | sequence |
---|---|---|---|---|---|
chr20 | 189774 | 189792 | + | 0.8794931 | ACTCCTCTAGAGGGTGCTC |
chr20 | 239773 | 239791 | + | 0.9003697 | TTGCCACTGGGGGGAGACA |
chr20 | 247783 | 247801 | - | 0.9337214 | CTGCCGGCAGATGGCGGTA |
chr20 | 281074 | 281092 | - | 0.8511201 | TTGCCTGCAGGGGTGGGAA |
The interaction between TF and DNA would leave a “footprint” in motif position, but it is not evident in a single site, so integrated footprint is necessary. In addition, we only consider Tn5 cut site. This function is based on the motif scan.
First, collecting all cut site from the bed file(Note: every line in the bed file is a DNA fragment) and save them.
## extract cut site position from bed file
fra_path <- system.file("extdata", "chr20.50000.bed.bz2", package="esATAC")
frag <- as.vector(bunzip2(filename = fra_path,
destname = file.path(getwd(), "chr20.50000.bed"),
ext="bz2", FUN=bzfile, overwrite=TRUE, remove = FALSE))
cs.data <- extractcutsite(bedInput = frag, prefix = "ATAC")
Next, plot footprint for different motifs.
In the motif scan, we get a variable named “motif.data”, is contains multi-motif information. In order to plot footprint of these motif in a single procedure, we will use the output of function motifscan, here is “motif.data”.
fp <- atacCutSiteCount(atacProcCutSite = cs.data, atacProcMotifScan = motif.data)
The following is CTCF footprint using example data.
Note: we only using a small part of the chromatin 20 as example.
esATAC is organized in data flow graph. Except for referring manual, the user may print the map to know workflow order.
bedbzfile <- system.file(package="esATAC", "extdata", "chr20.50000.bed.bz2")
bedfile <- file.path(td,"chr20.50000.bed")
bunzip2(bedbzfile,destname=bedfile,overwrite=TRUE,remove=FALSE)
peakproc <-bedUtils(bedInput = bedfile,maxFragLen = 100, chrFilterList = NULL) %>%
atacPeakCalling
peakproc %>% printMap
By printing the map, it is easy to know what valid processes are available to call next and what preprocess has been done before.
It is easy to query the parameters set for the process with ATACProc objects. You can query available parameters like this:
#query all of available parameters
getParamItems(peakproc)
## [1] "bedInput" "bedFileList" "inBedDir" "bedOutput"
## [5] "outTmpDir" "fragmentSize" "fileformat" "verbose"
The value of a specific parameter can be obtain like this:
#query a parameter value
getParam(peakproc,"fragmentSize")
## [1] 0
Similarly, it is also easy to query the report value calculated by the process with ATACProc objects
sambzfile <- system.file(package="esATAC", "extdata", "Example.sam.bz2")
samfile <- file.path(td,"Example.sam")
bunzip2(sambzfile,destname=samfile,overwrite=TRUE,remove=FALSE)
samToBedProc<-samToBed(samInput = samfile)
When the ATACProc objects are obtained, you can query all of available report items.
getReportItems(samToBedProc)
## [1] "report" "total"
## [3] "save" "filted"
## [5] "extlen" "unique"
## [7] "multimap" "non-mitochondrial"
## [9] "non-mitochondrial-multimap"
The value of specific report item can be get like this:
#query a parameter value
getReportVal(samToBedProc,"report")
## Item Value
## 1 total 209
## 2 save 13
## 3 filted 150
## 4 extlen 0
## 5 unique 0
## 6 multimap 15
If the user call a process function that was called last time and finished, the process function will not redo the process. So if users need to redo the process, they have to clear the cache like this:
clearProcCache(peakproc)
process(peakproc)
We would like to thank Huan Fang for package testing and valuable suggestions,
and Kui Hua for providing package testing on Macbook.
[1] Langmead, B., & Salzberg, S. L. (2012). Fast gapped-read alignment with Bowtie 2. Nature methods, 9(4), 357-359.
[2] Schubert, Lindgreen, and Orlando (2016). AdapterRemoval v2: rapid adapter trimming, identification, and read merging. BMC Research Notes, 12;9(1):88.
[3] Boyle, A. P., Guinney, J., Crawford, G. E., & Furey, T. S. (2008). F-Seq: a feature density estimator for high-throughput sequence tags. Bioinformatics, 24(21), 2537-2538.
[4] Yu G, Wang L and He Q (2015). “ChIPseeker: an R/Bioconductor package for ChIP peak annotation, comparison and visualization.” Bioinformatics, 31(14), pp. 2382-2383. doi: 10.1093/bioinformatics/btv145.
[5] Yu G, Wang L, Han Y and He Q (2012). “clusterProfiler: an R package for comparing biological themes among gene clusters.” OMICS: A Journal of Integrative Biology, 16(5), pp. 284-287. doi: 10.1089/omi.2011.0118.
[6] Pagès H, Aboyoun P, Gentleman R and DebRoy S (2017). Biostrings: String objects representing biological sequences, and matching algorithms. R package version 2.44.2.
[7] Buenrostro, J. D., Giresi, P. G., Zaba, L. C., Chang, H. Y., & Greenleaf, W. J. (2013). Transposition of native chromatin for fast and sensitive epigenomic profiling of open chromatin, DNA-binding proteins and nucleosome position. Nature methods, 10(12), 1213-1218.
sessionInfo()
## R version 3.4.2 (2017-09-28)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 16.04.3 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.6-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.6-bioc/R/lib/libRlapack.so
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_US.UTF-8 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=en_US.UTF-8
## [9] LC_ADDRESS=en_US.UTF-8 LC_TELEPHONE=en_US.UTF-8
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=en_US.UTF-8
##
## attached base packages:
## [1] stats4 parallel stats graphics grDevices utils datasets
## [8] methods base
##
## other attached packages:
## [1] clusterProfiler_3.6.0
## [2] DOSE_3.4.0
## [3] ChIPseeker_1.14.0
## [4] Rbowtie2_1.0.0
## [5] R.utils_2.5.0
## [6] R.oo_1.21.0
## [7] R.methodsS3_1.7.1
## [8] org.Hs.eg.db_3.4.2
## [9] TxDb.Hsapiens.UCSC.hg19.knownGene_3.2.2
## [10] GenomicFeatures_1.30.0
## [11] AnnotationDbi_1.40.0
## [12] Biobase_2.38.0
## [13] BSgenome.Hsapiens.UCSC.hg19_1.4.0
## [14] BSgenome_1.46.0
## [15] rtracklayer_1.38.0
## [16] magrittr_1.5
## [17] bindrcpp_0.2
## [18] esATAC_1.0.3
## [19] Rsamtools_1.30.0
## [20] Biostrings_2.46.0
## [21] XVector_0.18.0
## [22] QuasR_1.18.0
## [23] Rbowtie_1.18.0
## [24] GenomicRanges_1.30.0
## [25] GenomeInfoDb_1.14.0
## [26] IRanges_2.12.0
## [27] S4Vectors_0.16.0
## [28] BiocGenerics_0.24.0
##
## loaded via a namespace (and not attached):
## [1] backports_1.1.1 fastmatch_1.1-0
## [3] VGAM_1.0-4 plyr_1.8.4
## [5] igraph_1.1.2 lazyeval_0.2.1
## [7] splines_3.4.2 BiocParallel_1.12.0
## [9] ggplot2_2.2.1 gridBase_0.4-7
## [11] TFBSTools_1.16.0 digest_0.6.12
## [13] BiocInstaller_1.28.0 htmltools_0.3.6
## [15] GOSemSim_2.4.0 viridis_0.4.0
## [17] GO.db_3.4.2 gdata_2.18.0
## [19] memoise_1.1.0 JASPAR2016_1.6.0
## [21] readr_1.1.1 annotate_1.56.0
## [23] matrixStats_0.52.2 prettyunits_1.0.2
## [25] colorspace_1.3-2 blob_1.1.0
## [27] dplyr_0.7.4 RCurl_1.95-4.8
## [29] jsonlite_1.5 bindr_0.1
## [31] TFMPvalue_0.0.6 brew_1.0-6
## [33] VariantAnnotation_1.24.0 glue_1.2.0
## [35] gtable_0.2.0 zlibbioc_1.24.0
## [37] UpSetR_1.3.3 DelayedArray_0.4.1
## [39] Rook_1.1-1 scales_0.5.0
## [41] futile.options_1.0.0 DBI_0.7
## [43] Rcpp_0.12.13 plotrix_3.6-6
## [45] viridisLite_0.2.0 xtable_1.8-2
## [47] progress_1.1.2 bit_1.1-12
## [49] htmlwidgets_0.9 httr_1.3.1
## [51] DiagrammeR_0.9.2 fgsea_1.4.0
## [53] gplots_3.0.1 RColorBrewer_1.1-2
## [55] rJava_0.9-9 pkgconfig_2.0.1
## [57] XML_3.98-1.9 rlang_0.1.2
## [59] reshape2_1.4.2 munsell_0.4.3
## [61] tools_3.4.2 visNetwork_2.0.1
## [63] downloader_0.4 DirichletMultinomial_1.20.0
## [65] RSQLite_2.0 evaluate_0.10.1
## [67] stringr_1.2.0 yaml_2.1.14
## [69] knitr_1.17 bit64_0.9-7
## [71] caTools_1.17.1 purrr_0.2.4
## [73] KEGGREST_1.18.0 poweRlaw_0.70.1
## [75] DO.db_2.9 biomaRt_2.34.0
## [77] compiler_3.4.2 rstudioapi_0.7
## [79] rgexf_0.15.3 png_0.1-7
## [81] tibble_1.3.4 stringi_1.1.5
## [83] highr_0.6 futile.logger_1.4.3
## [85] GenomicFiles_1.14.0 lattice_0.20-35
## [87] CNEr_1.14.0 Matrix_1.2-11
## [89] data.table_1.10.4-3 bitops_1.0-6
## [91] qvalue_2.10.0 R6_2.2.2
## [93] latticeExtra_0.6-28 hwriter_1.3.2
## [95] RMySQL_0.10.13 ShortRead_1.36.0
## [97] KernSmooth_2.23-15 gridExtra_2.3
## [99] lambda.r_1.2 boot_1.3-20
## [101] gtools_3.5.0 assertthat_0.2.0
## [103] seqLogo_1.44.0 SummarizedExperiment_1.8.0
## [105] rprojroot_1.2 GenomicAlignments_1.14.0
## [107] GenomeInfoDbData_0.99.1 hms_0.3
## [109] influenceR_0.1.0 VennDiagram_1.6.17
## [111] grid_3.4.2 tidyr_0.7.2
## [113] rmarkdown_1.6 rvcheck_0.0.9