Figure 1

The Figure 1 is the showcase for geomeTriD package to present multiple genomic signals along with 3D models in different layout.

Load Libraries

library(geomeTriD)
library(geomeTriD.documentation)
library(GenomicRanges)
library(TxDb.Hsapiens.UCSC.hg19.knownGene)
library(trackViewer)

Prepare the annotation and genomic data

The data were downloaded from ENCODE for GM12878 with assembly hg19. The 3D structure model were downloaded from FLAMINGO results for GEO dataset with accession GSE63525.

chr <- 'chr4'
range <- as(seqinfo(TxDb.Hsapiens.UCSC.hg19.knownGene)[chr], 'GRanges')
prefix <- 'https://www.encodeproject.org/files/'
url <- c(#ctcf = 'ENCFF364OXN/@@download/ENCFF364OXN.bigWig',
         #smc3 = 'ENCFF235BXX/@@download/ENCFF235BXX.bigWig',
         #rad21 = 'ENCFF567EGK/@@download/ENCFF567EGK.bigWig',
         #H3K4me3 = 'ENCFF674QZB/@@download/ENCFF674QZB.bigWig',
         #H3K27me3 = 'ENCFF594HSG/@@download/ENCFF594HSG.bigWig',
         H3K9me3 = 'ENCFF776OVW/@@download/ENCFF776OVW.bigWig',
         #H3K4me1 = 'ENCFF682WPF/@@download/ENCFF682WPF.bigWig',
         H3K27ac = 'ENCFF180LKW/@@download/ENCFF180LKW.bigWig')
urls <- paste0(prefix, url)
names(urls) <- names(url)
colorCode <- c('active'='#EA262E', 'inactive'='#179281')
# download the files and import the signals
genomicSigs <- importSignalFromUrl(urls, range = range,
                                   cols = list(
                                     'H3K9me3'=c('#000011', 'blue'),
                                     'H3K27ac'=c('#111100', 'orange')),
                                   format = 'BigWig')
## adding rname 'https://www.encodeproject.org/files/ENCFF776OVW/@@download/ENCFF776OVW.bigWig'
## Error while performing HEAD request.
##    Proceeding without cache information.
## adding rname 'https://www.encodeproject.org/files/ENCFF180LKW/@@download/ENCFF180LKW.bigWig'
## Error while performing HEAD request.
##    Proceeding without cache information.
# import the A/B compartments
compartments <- rtracklayer::import('https://ftp.ncbi.nlm.nih.gov/geo/series/GSE63nnn/GSE63525/suppl/GSE63525%5FGM12878%5Fsubcompartments.bed.gz')
compartments$col <- ifelse(is.na(compartments$name), 'gray', 
                           ifelse(grepl('^A', compartments$name),
                                  colorCode['active'],
                                  colorCode['inactive']))
compartments$label <- compartments$name
compartments$type <- 'compartment'
# import the 3D model generated by FLAMINGO
chrs <- paste0('chr', c(1:22, 'X'))
gm12878 <- paste0('https://github.com/wangjr03/FLAMINGO/',
                  'raw/refs/heads/main/predictions/GM12878/',
                  chrs, '_5kb.txt')
names(gm12878) <- chrs
gm12878.gl <- importFLAMINGO(gm12878)

Plot the data by geomeTriD

# create a list of threeJsGeometry objects
objsC <- view3dStructure(gm12878.gl[[chr]], # 3D model by FLAMINGO
                         k=3, # 3D
                         feature.gr=compartments[seqnames(compartments) %in%
                                                   chr], # subset of the compartments
                         genomicSigs = genomicSigs, # the histon markers
                         reverseGenomicSigs = FALSE, lwd.maxGenomicSigs = 10,
                         show_coor = FALSE, label_gene = FALSE,#close the labels
                         renderer = 'none', # return a list of objects
                         resolution = 1)
## top vs bottom layout
for(i in seq_along(objsC)){
  if(grepl('backbone|compartment', names(objsC)[i])){
    objsC[[i]]$layer <- 'bottom' ## put compartment and backbone to bottom layer
  }
}
threeJsViewer(objsC, background = 'black')
#widget <- threeJsViewer(objsC, background = 'black', width = 1920, height = 1080)
#tempfile <- 'Fig1.part1.html'
#htmlwidgets::saveWidget(widget, file=tempfile)
#utils::browseURL(tempfile)


## side by side layout
for(i in seq_along(objsC)){
  if(grepl('backbone|compartment', names(objsC)[i])){
    objsC[[i]]$layer <- 'top' ## put compartment to top layer
  }
  if(grepl('H3K9me3|H3K27ac', names(objsC)[i])){
    objsC[[i]]$side <- 'right' ## put H3K9me3 and H3K27ac signals to right panel
  }
}
threeJsViewer(objsC,
              title = c('A/B compartment',
                        'H3K9me3 and H3K27Ac signals'),
              background = 'black')
# widget <- threeJsViewer(objsC,
#               title = c('A/B compartment',
#                         'H3K9me3 and H3K27Ac signals'),
#               background = 'black', width = 1920, height = 1080)
# tempfile <- 'Fig1.part2.html'
# htmlwidgets::saveWidget(widget, file=tempfile)
# utils::browseURL(tempfile)

SessionInfo

## R version 4.5.0 (2025-04-11)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.2 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=en_US.UTF-8    
##  [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] grid      stats4    stats     graphics  grDevices utils     datasets 
## [8] methods   base     
## 
## other attached packages:
##  [1] trackViewer_1.45.0                     
##  [2] TxDb.Hsapiens.UCSC.hg19.knownGene_3.2.2
##  [3] GenomicFeatures_1.61.3                 
##  [4] AnnotationDbi_1.71.0                   
##  [5] Biobase_2.69.0                         
##  [6] GenomicRanges_1.61.0                   
##  [7] GenomeInfoDb_1.45.4                    
##  [8] IRanges_2.43.0                         
##  [9] S4Vectors_0.47.0                       
## [10] BiocGenerics_0.55.0                    
## [11] generics_0.1.4                         
## [12] geomeTriD.documentation_0.0.3          
## [13] geomeTriD_1.3.11                       
## 
## loaded via a namespace (and not attached):
##   [1] strawr_0.0.92               RColorBrewer_1.1-3         
##   [3] rstudioapi_0.17.1           jsonlite_2.0.0             
##   [5] magrittr_2.0.3              farver_2.1.2               
##   [7] rmarkdown_2.29              fs_1.6.6                   
##   [9] BiocIO_1.19.0               ragg_1.4.0                 
##  [11] vctrs_0.6.5                 memoise_2.0.1              
##  [13] Rsamtools_2.25.0            RCurl_1.98-1.17            
##  [15] base64enc_0.1-3             htmltools_0.5.8.1          
##  [17] S4Arrays_1.9.1              progress_1.2.3             
##  [19] plotrix_3.8-4               curl_6.3.0                 
##  [21] Rhdf5lib_1.31.0             rhdf5_2.53.1               
##  [23] SparseArray_1.9.0           Formula_1.2-5              
##  [25] sass_0.4.10                 parallelly_1.45.0          
##  [27] bslib_0.9.0                 htmlwidgets_1.6.4          
##  [29] desc_1.4.3                  Gviz_1.53.0                
##  [31] httr2_1.1.2                 cachem_1.1.0               
##  [33] GenomicAlignments_1.45.0    igraph_2.1.4               
##  [35] lifecycle_1.0.4             pkgconfig_2.0.3            
##  [37] Matrix_1.7-3                R6_2.6.1                   
##  [39] fastmap_1.2.0               MatrixGenerics_1.21.0      
##  [41] future_1.58.0               aricode_1.0.3              
##  [43] clue_0.3-66                 digest_0.6.37              
##  [45] colorspace_2.1-1            textshaping_1.0.1          
##  [47] Hmisc_5.2-3                 RSQLite_2.4.1              
##  [49] filelock_1.0.3              progressr_0.15.1           
##  [51] httr_1.4.7                  abind_1.4-8                
##  [53] compiler_4.5.0              withr_3.0.2                
##  [55] bit64_4.6.0-1               backports_1.5.0            
##  [57] htmlTable_2.4.3             BiocParallel_1.43.3        
##  [59] DBI_1.2.3                   R.utils_2.13.0             
##  [61] biomaRt_2.65.0              MASS_7.3-65                
##  [63] rappdirs_0.3.3              DelayedArray_0.35.1        
##  [65] rjson_0.2.23                tools_4.5.0                
##  [67] foreign_0.8-90              future.apply_1.20.0        
##  [69] nnet_7.3-20                 R.oo_1.27.1                
##  [71] glue_1.8.0                  restfulr_0.0.15            
##  [73] dbscan_1.2.2                InteractionSet_1.37.0      
##  [75] rhdf5filters_1.21.0         checkmate_2.3.2            
##  [77] cluster_2.1.8.1             gtable_0.3.6               
##  [79] BSgenome_1.77.0             R.methodsS3_1.8.2          
##  [81] ensembldb_2.33.0            data.table_1.17.4          
##  [83] hms_1.1.3                   xml2_1.3.8                 
##  [85] XVector_0.49.0              RANN_2.6.2                 
##  [87] pillar_1.10.2               stringr_1.5.1              
##  [89] dplyr_1.1.4                 BiocFileCache_2.99.5       
##  [91] lattice_0.22-7              deldir_2.0-4               
##  [93] rtracklayer_1.69.0          bit_4.6.0                  
##  [95] biovizBase_1.57.0           tidyselect_1.2.1           
##  [97] Biostrings_2.77.1           knitr_1.50                 
##  [99] gridExtra_2.3               ProtGenerics_1.41.0        
## [101] SummarizedExperiment_1.39.0 xfun_0.52                  
## [103] matrixStats_1.5.0           stringi_1.8.7              
## [105] UCSC.utils_1.5.0            lazyeval_0.2.2             
## [107] yaml_2.3.10                 evaluate_1.0.3             
## [109] codetools_0.2-20            interp_1.1-6               
## [111] tibble_3.3.0                cli_3.6.5                  
## [113] rpart_4.1.24                systemfonts_1.2.3          
## [115] jquerylib_0.1.4             dichromat_2.0-0.1          
## [117] Rcpp_1.0.14                 globals_0.18.0             
## [119] grImport_0.9-7              dbplyr_2.5.0               
## [121] png_0.1-8                   XML_3.99-0.18              
## [123] parallel_4.5.0              pkgdown_2.1.3              
## [125] rgl_1.3.18                  ggplot2_3.5.2              
## [127] blob_1.2.4                  prettyunits_1.2.0          
## [129] jpeg_0.1-11                 latticeExtra_0.6-30        
## [131] AnnotationFilter_1.33.0     bitops_1.0-9               
## [133] txdbmaker_1.5.4             listenv_0.9.1              
## [135] VariantAnnotation_1.55.0    scales_1.4.0               
## [137] purrr_1.0.4                 crayon_1.5.3               
## [139] rlang_1.1.6                 KEGGREST_1.49.0