Single-molecule assays like NOMe-seq, and dSMF are superior to DNase-seq and ATAC-seq as they do not cleave DNA (see schematics below). Thus, they enable quantification of all three i.e., protein-free, Transcription Factor-bound and histone-complex-bound states. But a user-friendly tool to visualize and quantify such states is lacking. SMTrackR, a Bioconductor package visualizes protein-DNA binding states on individual sequenced DNA molecules. SMTrackR queries the single-molecule footprint database we built and hosted at Galaxy Server. It comprises of BigBed files generated from NOMe-seq and dSMF datasets. SMTrackR exploits UCSC REST API to query a BigBed file and plot footprint heatmap categorized in different binding states as well as report their occupancies.
SMTrackR::plotFootprints() will generate the heatmap
pdf file stored in
<target_dir>/plots/<label>.heatmap.pdf if
target_dir option is provided, otherwise the plot will be
saved in the tempdir, the fuction reports the path in both
the cases. A Gviz compatible code can also be generated by
running the generateGvizCodeforSMF function. See the code
chunk below.
SMTrackR::plotFootprints(organism = "mmusculus", model = "8cell", condition = "WT",
genome_assembly = "mm10", type = "SMF", chromosome = "chr5",
start = 113847750, end = 113847780, tr = "8cell",
label = "8cell", fp_cap = 50, remove_dup = TRUE, target_dir = "smf")## Using UCSC API to get JSON dump from BigBed Track ...
## adding rname 'https://api.genome.ucsc.edu/getData/track?hubUrl=https://usegalaxy.org/api/datasets/f9cad7b01a472135f6860d3104863246/display?to_ext=txt;genome=mm10;track=8cell;chrom=chr5;start=113847750;end=113847780'
## Heatmap is saved in the file smf/plots/8cell.heatmap.pdf
##
## [1] "smf/plots/8cell.heatmap.pdf"
SMTrackR::generateGvizCodeforSMF(organism = "mmusculus", model = "8cell", condition = "WT",
genome_assembly = "mm10", type = "SMF", chromosome = "chr5",
start = 113847750, end = 113847780, tr = "8cell",
label = "8cell", fp_cap = 50, remove_dup = TRUE,
gviz_left = 1000, gviz_right = 1000, target_dir = "smf")## A R code file smf/8cell.gviz.R is written. Please run the command: Rscript smf/8cell.gviz.R to place the heatmap on gene track!
To learn more about dSMF data, please read the article Krebs et al., Mol. Cell., 2017. Briefly, it uses two exogenous methyltransferases in tandem to map protein-DNA binding on individual sequenced DNA molecules.
The command below will generate plots/peak229.plot.pdf
file.
## Using UCSC API to get JSON dump from BigBed Track ...
## adding rname 'https://api.genome.ucsc.edu/getData/track?hubUrl=https://usegalaxy.org/api/datasets/f9cad7b01a472135e6897af4f6ed17e1/display?to_ext=txt;genome=dm6;track=fp_and_mvec;chrom=chr2L;start=480290;end=480320'
## Heatmap is saved in the file dsmf/plots/peak229.heatmap.pdf
##
## [1] "dsmf/plots/peak229.heatmap.pdf"
Running the command below will generate a pdf file
plots/smac_seq.plot.pdf.
## Using UCSC API to get JSON dump from BigBed Track ...
## adding rname 'https://api.genome.ucsc.edu/getData/track?hubUrl=https://usegalaxy.org/api/datasets/f9cad7b01a472135c1c27ee41f61cd34/display?to_ext=txt;genome=sacCer3;track=nanopore_meth_calls;chrom=chrIII;start=114300;end=114600'
## Heatmap is saved in the file nanopore/plots/smac_seq.heatmap.pdf
##
All uptodate tracks are available here and can be listed using the following command.
## R version 4.5.2 (2025-10-31)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.4 LTS
##
## Matrix products: default
## BLAS: /usr/lib/x86_64-linux-gnu/openblas-pthread/libblas.so.3
## LAPACK: /usr/lib/x86_64-linux-gnu/openblas-pthread/libopenblasp-r0.3.26.so; LAPACK version 3.12.0
##
## 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
##
## time zone: Etc/UTC
## tzcode source: system (glibc)
##
## attached base packages:
## [1] stats graphics grDevices utils datasets methods base
##
## other attached packages:
## [1] SMTrackR_0.99.6 BiocStyle_2.39.0
##
## loaded via a namespace (and not attached):
## [1] SummarizedExperiment_1.41.1 httr2_1.2.2
## [3] rjson_0.2.23 xfun_0.56
## [5] bslib_0.10.0 Biobase_2.71.0
## [7] lattice_0.22-9 vctrs_0.7.1
## [9] tools_4.5.2 bitops_1.0-9
## [11] generics_0.1.4 stats4_4.5.2
## [13] curl_7.0.0 parallel_4.5.2
## [15] tibble_3.3.1 RSQLite_2.4.6
## [17] blob_1.3.0 pkgconfig_2.0.3
## [19] Matrix_1.7-4 dbplyr_2.5.2
## [21] S4Vectors_0.49.0 cigarillo_1.1.0
## [23] lifecycle_1.0.5 compiler_4.5.2
## [25] stringr_1.6.0 Rsamtools_2.27.1
## [27] Biostrings_2.79.5 Seqinfo_1.1.0
## [29] codetools_0.2-20 htmltools_0.5.9
## [31] sys_3.4.3 buildtools_1.0.0
## [33] sass_0.4.10 RCurl_1.98-1.17
## [35] yaml_2.3.12 pillar_1.11.1
## [37] crayon_1.5.3 jquerylib_0.1.4
## [39] BiocParallel_1.45.0 cachem_1.1.0
## [41] DelayedArray_0.37.0 abind_1.4-8
## [43] tidyselect_1.2.1 digest_0.6.39
## [45] stringi_1.8.7 purrr_1.2.1
## [47] dplyr_1.2.0 restfulr_0.0.16
## [49] maketools_1.3.2 fastmap_1.2.0
## [51] grid_4.5.2 cli_3.6.5
## [53] SparseArray_1.11.11 magrittr_2.0.4
## [55] S4Arrays_1.11.1 XML_3.99-0.22
## [57] withr_3.0.2 filelock_1.0.3
## [59] rappdirs_0.3.4 bit64_4.6.0-1
## [61] rmarkdown_2.30 XVector_0.51.0
## [63] httr_1.4.8 matrixStats_1.5.0
## [65] bit_4.6.0 memoise_2.0.1
## [67] evaluate_1.0.5 knitr_1.51
## [69] GenomicRanges_1.63.1 IRanges_2.45.0
## [71] BiocIO_1.21.0 BiocFileCache_3.1.0
## [73] rtracklayer_1.71.3 rlang_1.1.7
## [75] glue_1.8.0 DBI_1.3.0
## [77] BiocManager_1.30.27 BiocGenerics_0.57.0
## [79] jsonlite_2.0.0 R6_2.6.1
## [81] MatrixGenerics_1.23.0 GenomicAlignments_1.47.0