Dandelion Plot

Sometimes, there are as many as hundreds of SNPs invoved in one gene. Dandelion plot can be used to depict such dense SNPs. Please note that the height of the dandelion indicates the desity of the SNPs.

library(trackViewer)
library(TxDb.Hsapiens.UCSC.hg19.knownGene)
library(org.Hs.eg.db)
library(rtracklayer)
methy <- import(system.file("extdata", "methy.bed", package="trackViewer"), "BED")
gr <- GRanges("chr22", IRanges(50968014, 50970514, names="TYMP"))
trs <- geneModelFromTxdb(TxDb.Hsapiens.UCSC.hg19.knownGene,
                         org.Hs.eg.db,
                         gr=gr)
features <- c(range(trs[[1]]$dat), range(trs[[5]]$dat))
names(features) <- c(trs[[1]]$name, trs[[5]]$name)
features$fill <- c("lightblue", "mistyrose")
features$height <- c(.02, .04)
dandelion.plot(methy, features, ranges=gr, type="pin")

Change the type of Dandelion plot

There are one more type for dandelion plot, i.e., type “fan”. The area of the fan indicates the percentage of methylation or rate of mutation.

methy$color <- 3
methy$border <- "gray"
## Score info is required and the score must be a number in [0, 1]
m <- max(methy$score)
methy$score <- methy$score/m
dandelion.plot(methy, features, ranges=gr, type="fan")

methy$color <- rep(list(c(3, 5)), length(methy))
methy$score2 <- (max(methy$score) - methy$score)/m
legends <- list(list(labels=c("s1", "s2"), fill=c(3, 5)))
dandelion.plot(methy, features, ranges=gr, type="pie", legend=legends)

Change the number of dandelions

## Less dandelions
dandelion.plot(methy, features, ranges=gr, type="circle", maxgaps=1/10)

## More dandelions
dandelion.plot(methy, features, ranges=gr, type="circle", maxgaps=1/100)

Users can also specify the maximum distance between neighboring dandelions by setting the maxgaps as a GRanges object.

maxgaps <- tile(gr, n = 10)[[1]]
dandelion.plot(methy, features, ranges=gr, type="circle", maxgaps=maxgaps)

Add y-axis (yaxis)

Set yaxis to TRUE to add y-axis, and set heightMethod=mean to use the mean score as the height.

dandelion.plot(methy, features, ranges=gr, type="pie", 
               maxgaps=1/100, yaxis = TRUE, heightMethod = mean,
               ylab='mean of methy scores')

yaxis = c(0, 0.5, 1)
dandelion.plot(methy, features, ranges=gr, type="pie", 
               maxgaps=1/100, yaxis = yaxis, heightMethod = mean,
               ylab='mean of methy scores')

Session Info

sessionInfo()

R version 4.5.0 Patched (2025-04-21 r88169) Platform: aarch64-apple-darwin20 Running under: macOS Ventura 13.7.1

Matrix products: default BLAS: /Library/Frameworks/R.framework/Versions/4.5-arm64/Resources/lib/libRblas.0.dylib LAPACK: /Library/Frameworks/R.framework/Versions/4.5-arm64/Resources/lib/libRlapack.dylib; LAPACK version 3.12.1

locale: [1] C/en_US.UTF-8/en_US.UTF-8/C/en_US.UTF-8/en_US.UTF-8

time zone: America/New_York tzcode source: internal

attached base packages: [1] grid stats4 stats graphics grDevices utils datasets [8] methods base

other attached packages: [1] httr_1.4.7
[2] VariantAnnotation_1.55.0
[3] Rsamtools_2.25.0
[4] Biostrings_2.77.0
[5] XVector_0.49.0
[6] SummarizedExperiment_1.39.0
[7] MatrixGenerics_1.21.0
[8] matrixStats_1.5.0
[9] org.Hs.eg.db_3.21.0
[10] TxDb.Hsapiens.UCSC.hg19.knownGene_3.2.2 [11] GenomicFeatures_1.61.0
[12] AnnotationDbi_1.71.0
[13] Biobase_2.69.0
[14] Gviz_1.53.0
[15] rtracklayer_1.69.0
[16] trackViewer_1.45.0
[17] GenomicRanges_1.61.0
[18] GenomeInfoDb_1.45.0
[19] IRanges_2.43.0
[20] S4Vectors_0.47.0
[21] BiocGenerics_0.55.0
[22] generics_0.1.3

loaded via a namespace (and not attached): [1] RColorBrewer_1.1-3 strawr_0.0.92 rstudioapi_0.17.1
[4] jsonlite_2.0.0 magrittr_2.0.3 farver_2.1.2
[7] rmarkdown_2.29 BiocIO_1.19.0 vctrs_0.6.5
[10] memoise_2.0.1 RCurl_1.98-1.17 base64enc_0.1-3
[13] htmltools_0.5.8.1 S4Arrays_1.9.0 progress_1.2.3
[16] curl_6.2.2 Rhdf5lib_1.31.0 SparseArray_1.9.0
[19] Formula_1.2-5 rhdf5_2.53.0 sass_0.4.10
[22] bslib_0.9.0 htmlwidgets_1.6.4 httr2_1.1.2
[25] cachem_1.1.0 GenomicAlignments_1.45.0 lifecycle_1.0.4
[28] pkgconfig_2.0.3 Matrix_1.7-3 R6_2.6.1
[31] fastmap_1.2.0 GenomeInfoDbData_1.2.14 digest_0.6.37
[34] colorspace_2.1-1 Hmisc_5.2-3 RSQLite_2.3.9
[37] filelock_1.0.3 abind_1.4-8 compiler_4.5.0
[40] bit64_4.6.0-1 htmlTable_2.4.3 backports_1.5.0
[43] BiocParallel_1.43.0 DBI_1.2.3 biomaRt_2.65.0
[46] rappdirs_0.3.3 DelayedArray_0.35.1 rjson_0.2.23
[49] tools_4.5.0 foreign_0.8-90 nnet_7.3-20
[52] glue_1.8.0 restfulr_0.0.15 InteractionSet_1.37.0
[55] rhdf5filters_1.21.0 checkmate_2.3.2 cluster_2.1.8.1
[58] gtable_0.3.6 BSgenome_1.77.0 ensembldb_2.33.0
[61] data.table_1.17.0 hms_1.1.3 xml2_1.3.8
[64] pillar_1.10.2 stringr_1.5.1 dplyr_1.1.4
[67] BiocFileCache_2.99.0 lattice_0.22-7 bit_4.6.0
[70] deldir_2.0-4 biovizBase_1.57.0 tidyselect_1.2.1
[73] knitr_1.50 gridExtra_2.3 ProtGenerics_1.41.0
[76] xfun_0.52 stringi_1.8.7 UCSC.utils_1.5.0
[79] lazyeval_0.2.2 yaml_2.3.10 evaluate_1.0.3
[82] codetools_0.2-20 interp_1.1-6 tibble_3.2.1
[85] BiocManager_1.30.25 cli_3.6.5 rpart_4.1.24
[88] jquerylib_0.1.4 dichromat_2.0-0.1 Rcpp_1.0.14
[91] grImport_0.9-7 dbplyr_2.5.0 png_0.1-8
[94] XML_3.99-0.18 parallel_4.5.0 ggplot2_3.5.2
[97] blob_1.2.4 prettyunits_1.2.0 latticeExtra_0.6-30
[100] jpeg_0.1-11 AnnotationFilter_1.33.0 bitops_1.0-9
[103] txdbmaker_1.5.0 scales_1.4.0 crayon_1.5.3
[106] BiocStyle_2.37.0 rlang_1.1.6 KEGGREST_1.49.0