ALPS (AnaLysis routines for ePigenomicS data) is an R package (Venu 2019) that provides tools for the analysis and to produce publication-ready visualizations, mainly aimed at genome-wide epigenomics data, e.g. ChIP-seq, ATAC-seq etc.
Bigwig files evolved to be a multi-purpose compressed binary format to store genome-wide data at base pair level. Bigwig files are mostly used to store genome-wide quantitative data such as ChIP-seq, ATAC-seq, WGBS, GRO-seq etc. Following figure illsutrates the important usecases with bigwig files.
There are multiple ways one can generate bigwig files from BAM files, using UCSC kent utils (ucscGenomeBrowser 2019) or with the deeptools bamCoverage function (Ramírez et al. 2014), which is the easiest way. Once the normalized bigwig files are generated and peaks are identified from BAM files, one would seldom use BAM files again in the entire workflow. The requirements of all downstream processes can be satisified with normalized bigwig files, e.g quantifying normalized read counts at peaks or promoters, visualizing enrichments in genome broswer or igv.
After the peaks are identified, the immediate steps would be to quantify normalized read counts at the identified peaks in order to perform explorative data analysis (EDA), PCA, unsupervised clustering to identify patterns among samples under consideration and generate novel biological insights.
ALPS
package is designed in a way to start with a minimal set of input and to reach a rich source of insights from the data. At the most, most functions in ALPS
require a data table with paths to bigwig files and associated sample meta information. Various functions will utilize this data table and generate downstream outputs. The package produces publication quality visualizations, of which most can be customized within R using ggplot2
ecosystem.
Following is the overview of the ALPS
workflow and available functions
Install the ALPS
package with the following code
To demonstrate the utility of different functions, ALPS
package comes with a set of example files that were taken from TCGA consortium’s ATAC-seq data from here published at (Corces et al. 2018).
Following steps walk you through loading the example data and how to use different function and how to integrate function’s output with other R/bioconductor packages to ease the workflow process.
## load the library
library(ALPS)
#>
#>
#> Registered S3 method overwritten by 'GGally':
#> method from
#> +.gg ggplot2
Most of the explorative analyses in epigenomics starts with quantifying enrichments or methylations at a set of genomic regions e.g. promoter regions or identified peak regions. This quantifications will be used as an input to downstream analyses such as PCA, clustering. The function multiBigwig_summary
takes sample data table with bigwig paths and corresponding bed file paths calculates enrichments at genomic regions. This function is a wrapper around rtracklayer
bigwig utilities (Lawrence, Gentleman, and Carey 2009). The function simultaneously generates the consensus peak-set from all the bed files present in input data table before calculating enrichments.
Read data table from ALPS
package
chr21_data_table <- system.file("extdata/bw", "ALPS_example_datatable.txt", package = "ALPS", mustWork = TRUE)
## attach path to bw_path and bed_path
d_path <- dirname(chr21_data_table)
chr21_data_table <- read.delim(chr21_data_table, header = TRUE)
chr21_data_table$bw_path <- paste0(d_path, "/", chr21_data_table$bw_path)
chr21_data_table$bed_path <- paste0(d_path, "/", chr21_data_table$bed_path)
chr21_data_table %>% head
#> sample_id group color_code
#> 1 ACCx_1 ACCx #8DD3C7
#> 2 ACCx_2 ACCx #8DD3C7
#> 3 BRCA_1 BRCA #FFED6F
#> 4 BRCA_2 BRCA #FFED6F
#> 5 GBMx_1 GBMx #B3DE69
#> 6 GBMx_2 GBMx #B3DE69
#> bw_path
#> 1 /tmp/RtmpeKrirc/Rinst17a47c6458d/ALPS/extdata/bw/ACCx_025FE5F8_T1.bw
#> 2 /tmp/RtmpeKrirc/Rinst17a47c6458d/ALPS/extdata/bw/ACCx_025FE5F8_T2.bw
#> 3 /tmp/RtmpeKrirc/Rinst17a47c6458d/ALPS/extdata/bw/BRCA_000CFD9F_T2.bw
#> 4 /tmp/RtmpeKrirc/Rinst17a47c6458d/ALPS/extdata/bw/BRCA_01112370_T1.bw
#> 5 /tmp/RtmpeKrirc/Rinst17a47c6458d/ALPS/extdata/bw/GBMx_09C0DCE7_T1.bw
#> 6 /tmp/RtmpeKrirc/Rinst17a47c6458d/ALPS/extdata/bw/GBMx_09C0DCE7_T2.bw
#> bed_path
#> 1 /tmp/RtmpeKrirc/Rinst17a47c6458d/ALPS/extdata/bw/ACCx_1.bed
#> 2 /tmp/RtmpeKrirc/Rinst17a47c6458d/ALPS/extdata/bw/ACCx_2.bed
#> 3 /tmp/RtmpeKrirc/Rinst17a47c6458d/ALPS/extdata/bw/BRCA_1.bed
#> 4 /tmp/RtmpeKrirc/Rinst17a47c6458d/ALPS/extdata/bw/BRCA_2.bed
#> 5 /tmp/RtmpeKrirc/Rinst17a47c6458d/ALPS/extdata/bw/GBMx_1.bed
#> 6 /tmp/RtmpeKrirc/Rinst17a47c6458d/ALPS/extdata/bw/GBMx_2.bed
Now run the function multiBigwig_summary
to calculate enrichments from all bigwig files by simultaneosly preparing consensus peak-set from all bed files in the column bed_path
enrichments <- multiBigwig_summary(data_table = chr21_data_table,
summary_type = "mean",
parallel = FALSE)
enrichments %>% head
#> chr start end ACCx_1 ACCx_2 BRCA_1 BRCA_2
#> 1 chr21 45639550 45640549 3.4293169 4.30173375 6.229088 6.053632
#> 2 chr21 45640994 45641493 5.1058709 4.39355741 2.412814 3.262252
#> 3 chr21 45642859 45643914 118.9869371 132.23839670 63.987536 82.130318
#> 4 chr21 45668163 45668662 0.3810360 0.02440864 4.359150 1.732686
#> 5 chr21 45689906 45690405 0.4229491 0.81362319 2.638012 0.923278
#> 6 chr21 45698386 45698910 46.5824657 54.70859597 1.141300 1.746901
#> GBMx_1 GBMx_2 LGGx_1 LGGx_2
#> 1 6.2633310 8.964782 8.0253673 6.3168158
#> 2 6.2341950 5.569840 2.2929616 4.2814719
#> 3 23.8171352 29.146589 57.2475263 50.2498071
#> 4 0.6438128 1.015109 0.8307302 0.9188698
#> 5 1.1011821 0.723354 2.4658742 1.0604877
#> 6 24.9284784 21.264368 2.8639639 2.1893534
With little teaking, the output from multiBigwig_summary
can be very easily integrated with other R/bioconductor packages for explorative analysis, PCA or clustering.
The function get_variable_regions
takes the output of multiBigwig_summary
or a similar format and returns a n
number of scaled variable regions, which can directly be used with tools such as ComplexHeatmap
(Gu, Eils, and Schlesner 2016).
Following is an example on how to integrate multiBigwig_summary
output to ComplexHeatmap
via get_variable_regions
enrichments_matrix <- get_variable_regions(enrichments_df = enrichments,
log_transform = TRUE,
scale = TRUE,
num_regions = 100)
suppressPackageStartupMessages(require(ComplexHeatmap))
suppressPackageStartupMessages(require(circlize))
Heatmap(enrichments_matrix, name = "enrichments",
col = colorRamp2(c(-1, 0, 1), c("green", "white", "red")),
show_row_names = FALSE,
show_column_names = TRUE,
show_row_dend = FALSE,
column_names_gp = gpar(fontsize = 8))
It is often of high interest in genome-wide quantitative data to check the correlations among replicates within a subgroup to identify specific patterns in the data. plot_correlation
function is designed for such use cases. The function is compatible with the output of multiBigwig_summary
and also with other tools output with similar format
plot_correlation(enrichments_df = enrichments,
log_transform = TRUE,
plot_type = "replicate_level",
sample_metadata = chr21_data_table)
Instead of correlations of replicates within and across groups, one can also do group level correlations after averaging all samples within a group. The argument plot_type = "group_level"
in plot_correlation
exactly does this.
## group_level
plot_correlation(enrichments_df = enrichments,
log_transform = TRUE,
plot_type = "group_level",
sample_metadata = chr21_data_table)
Either replicate_level
or group_level
plot appearance can be further modified with arguments that passed to corrplot::corrplot
or GGally::ggpairs
respectively.
Once the group-specific genomic regions (or peaks) identified with various differential enrichments packages, e.g. DESeq2, diffBind, QSEA, one would be interested to visualize enrichment qunatities across all samples of all groups to show magnitude of differnce in enrichments. plot_enrichments
function takes a data.frame of enrichments, either the output from multiBigwig_summary
or a similar format and plots enrichments in a combination of box and violin plots. The function is motivated by the paper (Allen M 2018) and a ggplot2 extension gghalves
(Frederik 2019). There are two ways one can plot enrichment differences, one way is to directly plot group level enrichments after averaging all samples within a group for each region and the other way is plotting paired conditions for each group, e.g. untreated, treated enrichments for a transcription factor. In both cases function needs a sample_metadata
table along with the enrichments
data.frame.
Following example illustrates the uses of plot_enrichments
function uses in different settings. If plot_type = "separate"
, function plots group level enrichments
## plot_type == "separate"
plot_enrichments(enrichments_df = enrichments,
log_transform = TRUE,
plot_type = "separate",
sample_metadata = chr21_data_table)
If plot_type = "overlap"
, function plots box plots along with overlap violins to show the distributions in paired conditions. The sample_metadata
for these plots require one more additional column which describes sample status. See the following example
## plot_type == "overlap"
enrichemnts_4_overlapviolins <- system.file("extdata/overlap_violins",
"enrichemnts_4_overlapviolins.txt",
package = "ALPS", mustWork = TRUE)
enrichemnts_4_overlapviolins <- read.delim(enrichemnts_4_overlapviolins, header = TRUE)
## metadata associated with above enrichments
data_table_4_overlapviolins <- system.file("extdata/overlap_violins",
"data_table_4_overlapviolins.txt",
package = "ALPS", mustWork = TRUE)
data_table_4_overlapviolins <- read.delim(data_table_4_overlapviolins, header = TRUE)
## enrichments table
enrichemnts_4_overlapviolins %>% head
#> chr start end s1 s2 s3 s4 s5
#> 1 chr21 5101659 5102227 0.25526578 0.4611994 0.21438043 0.1573575 -0.05081961
#> 2 chr21 5128223 5128741 0.82057469 0.9359073 0.66987324 0.8718381 0.31443426
#> 3 chr21 5154593 5155112 0.80378039 0.8452937 0.61775645 0.5620449 0.29282731
#> 4 chr21 5220488 5221005 0.04583463 0.7161867 0.04999978 0.5506898 0.41634277
#> 5 chr21 5221136 5221912 0.87553172 0.1238063 0.94939878 0.5923653 0.24742289
#> 6 chr21 5223327 5223826 0.12800909 0.5948275 0.58255203 0.7278129 -0.14432328
#> s6 s7 s8
#> 1 0.5747147 0.56832485 0.3590694
#> 2 1.1137688 -0.20168400 0.8110701
#> 3 0.8309985 -0.04550986 0.9127605
#> 4 0.8684639 0.37864897 0.3802162
#> 5 0.3206298 -0.13936118 0.7587155
#> 6 0.3536091 0.48731594 0.8615313
## metadata table
data_table_4_overlapviolins %>% head
#> bw_path sample_id sample_status group color_code
#> 1 path.bw s1 untreated tf1 gray
#> 2 path.bw s2 untreated tf2 gray
#> 3 path.bw s3 untreated tf1 gray
#> 4 path.bw s4 untreated tf2 gray
#> 5 path.bw s5 treated tf1 red
#> 6 path.bw s6 treated tf2 red
plot_enrichments(enrichments_df = enrichemnts_4_overlapviolins,
log_transform = FALSE,
plot_type = "overlap",
sample_metadata = data_table_4_overlapviolins,
overlap_order = c("untreated", "treated"))
There are additional arguments available for both separate
and overlap
to modify the appearance (please check ?plot_enrichments
), moreover the function returns a ggplot2
object which enables the user to change additional components of the plot.
In any genome-wide epigenomic analyses, it is often interesting to check enrichments at certain genomic loci, e.g. various histone modifications at a genomic region that define a chromatin state or co-binding of different transcription factors at a promoter or enhancer element. The classical way to acheive this task is to load all bigwig files into IGV or create a data hub at UCSC and navigate to the region of interest. This is not always practical and needs a substantial manual effort, in addition one requires a UCSC genome browser server in order to get this task done with the unpublished data.
To circumvent this problem, several R/bioconductor packages were designed (e.g. Gviz
, karyoploteR
). Even within R environment, one needs to put a significant effort to create UCSC genome browser like figures. The function plot_broswer_tracks
in ALPS
package requires a minimal input of a data table and a genomic region and produces a publication quality browser like plot. The function uses utilities within Gviz
package to generate the visualizations (Hahne and Ivanek 2016).
Following code snippet illustrates how one can use this function
One of the usual tasks in genome-wide epigenomic analyses is to identify the genomic locations of peaks/binding sites. This gives an overview of where a particular transcription factor frequently binds or where a particular type of histone modifications are observed. The function get_genomic_annotations
utilizes the above data table and returns the percentage of peaks or binding sites found in each of the genomic features such as promoters, UTRs, intergenetic regions etc. This function is a wrapper around ChIPseeker
’s annotatePeak
function (Yu, Wang, and He 2015).
Function also offers an option with merge_level
to merge overlaping peaks from different samples at different levels.
all
creates a consensus peak set by merging overlaping peaks from all samples present in the data_table
group_level
creates a group level consensus peak set. Meaning overlaping peaks from all samples of each group will be mergednone
does not create any consensus peak set. Per-sample genomic annotations will be returned
g_annotations <- get_genomic_annotations(data_table = chr21_data_table,
ref_gen = "hg38",
tss_region = c(-1000, 1000),
merge_level = "group_level")
#> Warning in data.table::dcast(., Feature ~ .id, value.var = "Frequency"): The
#> dcast generic in data.table has been passed a data.frame and will attempt to
#> redirect to the reshape2::dcast; please note that reshape2 is deprecated, and
#> this redirection is now deprecated as well. Please do this redirection yourself
#> like reshape2::dcast(.). In the next version, this warning will become an error.
g_annotations %>% head
#> Feature ACCx BRCA GBMx LGGx
#> 1 Promoter 45.762712 29.7709924 39.58333 40.000000
#> 2 5' UTR 1.694915 0.7633588 NA NA
#> 3 3' UTR 6.779661 3.8167939 6.25000 5.000000
#> 4 1st Exon 1.694915 2.2900763 6.25000 6.666667
#> 5 Other Exon 10.169492 9.9236641 6.25000 6.666667
#> 6 1st Intron 5.084746 3.0534351 NA 3.333333
The results returned from get_genomic_annotations
can directly be passed to the function plot_genomic_annotations
to visualize the percentages of peaks in each feature. The function can produce visualizations either as bar plot or heatmap
Transcription factors bind to DNA sequences with particular nucleotide stretches. A collection of all binding sites for a transcription factor can be represented by a sequence motif. Typically, transcription factor enrichment analyses utilizes these motif information and provide whether a given transcription factor is enriched in a given set of genomic regions, e.g. enhancers. Currently, there are number of different databases which provide transcription factor motifs in very different format. The function plot_motif_logo
takes input from various databases e.g. MEME, TRANSFAC, JASPAR, HOMER or a simple PFM and plots binding site representations either as a sequence logo or a barplot. Barplot representation has a several advantages over sequence logo which are described here. The logo plot utilizes the function ggseqlogo
from ggseqlogo
package (Wagih 2017).
Following example illustrates how to use the function to plot motif representations
ALPS
package benefited suggestions from
sessionInfo()
#> R version 4.0.3 (2020-10-10)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 18.04.5 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.12-bioc/R/lib/libRblas.so
#> LAPACK: /home/biocbuild/bbs-3.12-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=C
#> [9] LC_ADDRESS=C LC_TELEPHONE=C
#> [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
#>
#> attached base packages:
#> [1] grid stats graphics grDevices utils datasets methods
#> [8] base
#>
#> other attached packages:
#> [1] circlize_0.4.10 ComplexHeatmap_2.6.0 ALPS_1.4.0
#>
#> loaded via a namespace (and not attached):
#> [1] backports_1.1.10
#> [2] shadowtext_0.0.7
#> [3] ChIPseeker_1.26.0
#> [4] Hmisc_4.4-1
#> [5] fastmatch_1.1-0
#> [6] corrplot_0.84
#> [7] BiocFileCache_1.14.0
#> [8] plyr_1.8.6
#> [9] igraph_1.2.6
#> [10] lazyeval_0.2.2
#> [11] splines_4.0.3
#> [12] BiocParallel_1.24.0
#> [13] GenomeInfoDb_1.26.0
#> [14] ggplot2_3.3.2
#> [15] digest_0.6.27
#> [16] ensembldb_2.14.0
#> [17] htmltools_0.5.0
#> [18] GOSemSim_2.16.0
#> [19] magick_2.5.0
#> [20] viridis_0.5.1
#> [21] GO.db_3.12.0
#> [22] checkmate_2.0.0
#> [23] magrittr_1.5
#> [24] memoise_1.1.0
#> [25] BSgenome_1.58.0
#> [26] cluster_2.1.0
#> [27] annotate_1.68.0
#> [28] Biostrings_2.58.0
#> [29] graphlayouts_0.7.1
#> [30] matrixStats_0.57.0
#> [31] askpass_1.1
#> [32] enrichplot_1.10.0
#> [33] prettyunits_1.1.1
#> [34] jpeg_0.1-8.1
#> [35] colorspace_1.4-1
#> [36] blob_1.2.1
#> [37] rappdirs_0.3.1
#> [38] ggrepel_0.8.2
#> [39] xfun_0.18
#> [40] dplyr_1.0.2
#> [41] crayon_1.3.4
#> [42] RCurl_1.98-1.2
#> [43] TxDb.Hsapiens.UCSC.hg19.knownGene_3.2.2
#> [44] scatterpie_0.1.5
#> [45] genefilter_1.72.0
#> [46] VariantAnnotation_1.36.0
#> [47] survival_3.2-7
#> [48] glue_1.4.2
#> [49] polyclip_1.10-0
#> [50] gtable_0.3.0
#> [51] zlibbioc_1.36.0
#> [52] XVector_0.30.0
#> [53] GetoptLong_1.0.4
#> [54] TxDb.Hsapiens.UCSC.hg38.knownGene_3.10.0
#> [55] DelayedArray_0.16.0
#> [56] shape_1.4.5
#> [57] BiocGenerics_0.36.0
#> [58] scales_1.1.1
#> [59] DOSE_3.16.0
#> [60] DBI_1.1.0
#> [61] GGally_2.0.0
#> [62] Rcpp_1.0.5
#> [63] plotrix_3.7-8
#> [64] xtable_1.8-4
#> [65] htmlTable_2.1.0
#> [66] viridisLite_0.3.0
#> [67] progress_1.2.2
#> [68] clue_0.3-57
#> [69] foreign_0.8-80
#> [70] bit_4.0.4
#> [71] Formula_1.2-4
#> [72] stats4_4.0.3
#> [73] htmlwidgets_1.5.2
#> [74] httr_1.4.2
#> [75] fgsea_1.16.0
#> [76] gplots_3.1.0
#> [77] RColorBrewer_1.1-2
#> [78] ellipsis_0.3.1
#> [79] pkgconfig_2.0.3
#> [80] reshape_0.8.8
#> [81] XML_3.99-0.5
#> [82] farver_2.0.3
#> [83] ggseqlogo_0.1
#> [84] nnet_7.3-14
#> [85] Gviz_1.34.0
#> [86] dbplyr_1.4.4
#> [87] labeling_0.4.2
#> [88] tidyselect_1.1.0
#> [89] rlang_0.4.8
#> [90] reshape2_1.4.4
#> [91] AnnotationDbi_1.52.0
#> [92] munsell_0.5.0
#> [93] tools_4.0.3
#> [94] generics_0.0.2
#> [95] RSQLite_2.2.1
#> [96] evaluate_0.14
#> [97] stringr_1.4.0
#> [98] yaml_2.2.1
#> [99] org.Hs.eg.db_3.12.0
#> [100] knitr_1.30
#> [101] bit64_4.0.5
#> [102] tidygraph_1.2.0
#> [103] caTools_1.18.0
#> [104] purrr_0.3.4
#> [105] AnnotationFilter_1.14.0
#> [106] ggraph_2.0.3
#> [107] DO.db_2.9
#> [108] xml2_1.3.2
#> [109] biomaRt_2.46.0
#> [110] rstudioapi_0.11
#> [111] compiler_4.0.3
#> [112] curl_4.3
#> [113] png_0.1-7
#> [114] tibble_3.0.4
#> [115] tweenr_1.0.1
#> [116] stringi_1.5.3
#> [117] GenomicFeatures_1.42.0
#> [118] lattice_0.20-41
#> [119] ProtGenerics_1.22.0
#> [120] Matrix_1.2-18
#> [121] vctrs_0.3.4
#> [122] pillar_1.4.6
#> [123] lifecycle_0.2.0
#> [124] BiocManager_1.30.10
#> [125] GlobalOptions_0.1.2
#> [126] data.table_1.13.2
#> [127] cowplot_1.1.0
#> [128] bitops_1.0-6
#> [129] rtracklayer_1.50.0
#> [130] GenomicRanges_1.42.0
#> [131] qvalue_2.22.0
#> [132] R6_2.4.1
#> [133] latticeExtra_0.6-29
#> [134] KernSmooth_2.23-17
#> [135] gridExtra_2.3
#> [136] IRanges_2.24.0
#> [137] dichromat_2.0-0
#> [138] boot_1.3-25
#> [139] MASS_7.3-53
#> [140] gtools_3.8.2
#> [141] assertthat_0.2.1
#> [142] SummarizedExperiment_1.20.0
#> [143] rjson_0.2.20
#> [144] openssl_1.4.3
#> [145] GenomicAlignments_1.26.0
#> [146] Rsamtools_2.6.0
#> [147] S4Vectors_0.28.0
#> [148] GenomeInfoDbData_1.2.4
#> [149] parallel_4.0.3
#> [150] hms_0.5.3
#> [151] gghalves_0.1.0
#> [152] rpart_4.1-15
#> [153] tidyr_1.1.2
#> [154] rmarkdown_2.5
#> [155] rvcheck_0.1.8
#> [156] MatrixGenerics_1.2.0
#> [157] Cairo_1.5-12.2
#> [158] biovizBase_1.38.0
#> [159] ggforce_0.3.2
#> [160] Biobase_2.50.0
#> [161] base64enc_0.1-3
Allen M, Whitaker K, Poggiali D. 2018. “Raincloud plots: a multi-platform tool for robust data visualization.” PeerJ Preprints 6:E27137v1.
Corces, M. Ryan, Jeffrey M. Granja, Shadi Shams, Bryan H. Louie, Jose A. Seoane, Wanding Zhou, Tiago C. Silva, et al. 2018. “The Chromatin Accessibility Landscape of Primary Human Cancers.” Edited by Rehan Akbani, Christopher C. Benz, Evan A. Boyle, Bradley M. Broom, Andrew D. Cherniack, Brian Craft, John A. Demchok, et al. Science 362 (6413). American Association for the Advancement of Science. https://doi.org/10.1126/science.aav1898.
Frederik, Tiedemann. 2019. “gghalves: Easy half-half geoms in ggplot2.” https://cran.r-project.org/web/packages/gghalves.
Gu, Zuguang, Roland Eils, and Matthias Schlesner. 2016. “Complex heatmaps reveal patterns and correlations in multidimensional genomic data.” Bioinformatics 32 (18):2847–9. https://doi.org/10.1093/bioinformatics/btw313.
Hahne, Florian, and Robert Ivanek. 2016. “Visualizing Genomic Data Using Gviz and Bioconductor.” In Statistical Genomics: Methods and Protocols, edited by Sean Mathé Ewyand Davis, 335–51. New York, NY: Springer New York. https://doi.org/10.1007/978-1-4939-3578-9_16.
Lawrence, Michael, Robert Gentleman, and Vincent Carey. 2009. “rtracklayer: an R package for interfacing with genome browsers.” Bioinformatics 25 (14):1841–2. https://doi.org/10.1093/bioinformatics/btp328.
Ramírez, Fidel, Friederike Dündar, Sarah Diehl, Björn A. Grüning, and Thomas Manke. 2014. “deepTools: a flexible platform for exploring deep-sequencing data.” Nucleic Acids Research 42 (W1):W187–W191. https://doi.org/10.1093/nar/gku365.
ucscGenomeBrowser. 2019. “UCSC Genome Browser source tree.” https://github.com/ucscGenomeBrowser/kent.
Venu, Thatikonda. 2019. “ALPS: AnaLysis routines for ePigenomicS data.” https://github.com/itsvenu.
Wagih, Omar. 2017. “ggseqlogo: a versatile R package for drawing sequence logos.” Bioinformatics 33 (22):3645–7. https://doi.org/10.1093/bioinformatics/btx469.
Yu, Guangchuang, Li-Gen Wang, and Qing-Yu He. 2015. “ChIPseeker: An R/Bioconductor Package for Chip Peak Annotation, Comparison and Visualization.” Bioinformatics 31 (14):2382–3. https://doi.org/10.1093/bioinformatics/btv145.