1 Introduction

This document provides an introduction of the ELMER.data, which contains supporting data for ELMER (Yao, L., Shen, H., Laird, P. W., Farnham, P. J., & Berman, B. P. 2015). ELMER is package using DNA methylation to identify enhancers, and correlates enhancer state with expression of nearby genes to identify one or more transcriptional targets. Transcription factor (TF) binding site analysis of enhancers is coupled with expression analysis of all TFs to infer upstream regulators. ELMER.data provide 3 necessary data for ELMER analysis:

  1. Probes information: files with DNA methylation platforms metadata retrieved from http://zwdzwd.github.io/InfiniumAnnotation (Zhou, Wanding and Laird, Peter W and Shen, Hui 2016).
  2. Probes.motif: motif occurences within \(\pm 250bp\) of probe sites on HM450K/EPIC array aligned against hg19/hg38.

1.1 Installing and loading ELMER.data

To install this package, start R and enter

devtools::install_github(repo = "tiagochst/ELMER.data")
library("ELMER.data")
library("GenomicRanges")

2 Contents

2.1 Probes information

Probes information were retrieved from http://zwdzwd.github.io/InfiniumAnnotation (Zhou, Wanding and Laird, Peter W and Shen, Hui 2016).

for(plat in c("450K","EPIC")) {
  for(genome in c("hg38","hg19")) {
    base <- "http://zwdzwd.io/InfiniumAnnotation/current/"
    if(plat == "EPIC") {
      annotation <- paste0(base,"EPIC/EPIC.manifest.rda")
    } else {
      annotation <- paste0(base,"hm450/hm450.manifest.rda")
    }
    if(genome == "hg38") annotation <- gsub(".rda",".hg38.rda", annotation)
    if(!file.exists(basename(annotation))) {
      if(Sys.info()["sysname"] == "Windows") mode <- "wb" else  mode <- "w"
      downloader::download(annotation, basename(annotation), mode = mode)
    }
  }
}
data("EPIC.manifest")
as.data.frame(EPIC.manifest)[1:5,] %>% datatable(options = list(scrollX = TRUE,pageLength = 5)) 
data("EPIC.manifest.hg38")
as.data.frame(EPIC.manifest.hg38)[1:5,] %>% datatable(options = list(scrollX = TRUE,pageLength = 5)) 
data("hm450.manifest.hg38")
as.data.frame(hm450.manifest.hg38)[1:5,] %>% datatable(options = list(scrollX = TRUE,pageLength = 5)) 
data("hm450.manifest")
as.data.frame(hm450.manifest)[1:5,] %>% datatable(options = list(scrollX = TRUE,pageLength = 5)) 

2.2 TF family and subfamily classifications

ELMER uses the TFClass (Wingender, E., Schoeps, T., Haubrock, M., & Dönitz, J 2014), a classification of eukaryotic transcription factors based on the characteristics of their DNA-binding domains, to identify which are the TF that might be binding to the same region. For example, if a FOXA1 motif is found in a region, there is FOXA2 would also be able to bind to that region. For that ELMER uses two classifications, Family and sub-family. TFClass schema is shown below.

TFClass schema is below:

Level Rank denomination Definition Example
1 Superclass General topology of the DNA-binding domain Zinc-coordinating DNA-binding domains (Superclass 2)
2 Class Structural blueprint of the DNA-binding domain Nuclear receptors with C4 zinc fingers (Class 2.1)
3 Family Sequence & functional similarity Thyroid hormone receptor-related factors (NR1) (Family 2.1.2)
4 Subfamily Sequence-based subgroupings Retinoic acid receptors (NR1B) (Subfamily 2.1.2.1)
5 Genus Transcription factor gene RAR-α (Genus 2.1.2.1.1)
4 Species TF polypeptide RAR-α1 (Species 2.1.2.1.1.1)

The following code was used to create the objects:

library(xml2)
library(httr)
library(dplyr)
library(rvest)
createMotifRelevantTfs <- function(classification = "family"){
  message("Accessing hocomoco to get last version of TFs ", classification)
  file <- paste0(classification,".motif.relevant.TFs.rda")
   
    # Download from http://hocomoco.autosome.ru/human/mono
    tf.family <- "http://hocomoco11.autosome.ru/human/mono?full=true" %>% read_html()  %>%  html_table()
    tf.family <- tf.family[[1]]
    # Split TF for each family, this will help us map for each motif which are the some ones in the family
    # basicaly: for a TF get its family then get all TF in that family
    col <- ifelse(classification == "family", "TF family","TF subfamily")
    family <- split(tf.family,f = tf.family[[col]])
    
    motif.relevant.TFs <- plyr::alply(tf.family,1, function(x){  
      f <- x[[col]]
      if(f == "") return(x$`Transcription factor`) # Case without family, we will get only the same object
      return(unique(family[as.character(f)][[1]]$`Transcription factor`))
    },.progress = "text")
    #names(motif.relevant.TFs) <- tf.family$`Transcription factor`
    names(motif.relevant.TFs) <- tf.family$Model
    # Cleaning object
    attr(motif.relevant.TFs,which="split_type") <- NULL
    attr(motif.relevant.TFs,which="split_labels") <- NULL
  return(motif.relevant.TFs)
}
TF.family <- createMotifRelevantTfs("family")
TF.subfamily <- createMotifRelevantTfs("subfamily")
save(TF.family,file = "TF.family.rda", compress = "xz")
save(TF.subfamily,file = "TF.subfamily.rda", compress = "xz")
hocomoco.table <- "http://hocomoco11.autosome.ru/human/mono?full=true" %>% read_html()  %>%  html_table()
hocomoco.table <- hocomoco.table[[1]]
save(hocomoco.table,file = "data/hocomoco.table.rda", compress = "xz")

2.3 Probes.motif

Probes.motif provides information for motif occurences within\(\pm 250bp\) of probe sites on HM450K/EPIC array. HOMER (Heinz, Sven and Benner, Christopher and Spann, Nathanael and Bertolino, Eric and Lin, Yin C and Laslo, Peter and Cheng, Jason X and Murre, Cornelis and Singh, Harinder and Glass, Christopher K 2010) was used with a p-value < 0.0001 to scan a \(\pm 250bp\) region around each probe on HM450K/EPIC using HOCOMOCO V11 motif position weight matrices (PWMs) which provides transcription factor (TF) binding models for more than 600 human TFs. This data set is used in get.enriched.motif function in ELMER to calculate Odds Ratio of motif enrichments for a given set of probes. This data is storaged in a sparse matrix with wuth 640 columns, there is one matrix for HM450K aligned to hg19, one for HM450K aligned to hg38, one for EPIC aligned to hg19, one for EPIC aligned to hg38. Each row is each probe regions (annotation of the regions used can be found in this repository) and each column is motif from http://hocomoco.autosome.ru/ (Kulakovskiy, Ivan V and Medvedeva, Yulia A and Schaefer, Ulf and Kasianov, Artem S and Vorontsov, Ilya E and Bajic, Vladimir B and Makeev, Vsevolod 2013). The value 1 indicates the occurrence of a motif in a particular probe and 0 means no occurrence.

data("Probes.motif.hg19.450K")
dim(Probes.motif.hg19.450K)
## Loading required package: Matrix
## 
## Attaching package: 'Matrix'
## The following object is masked from 'package:S4Vectors':
## 
##     expand
## [1] 435391    771
str(Probes.motif.hg19.450K)
## Formal class 'ngCMatrix' [package "Matrix"] with 5 slots
##   ..@ i       : int [1:51917377] 1 2 11 18 20 23 27 29 35 36 ...
##   ..@ p       : int [1:772] 0 161636 175786 191191 206256 222389 265700 308214 348205 434144 ...
##   ..@ Dim     : int [1:2] 435391 771
##   ..@ Dimnames:List of 2
##   .. ..$ : chr [1:435391] "cg10029259" "cg14851608" "cg22163199" "cg09372808" ...
##   .. ..$ : chr [1:771] "AHR_HUMAN.H11MO.0.B" "AIRE_HUMAN.H11MO.0.C" "ALX1_HUMAN.H11MO.0.B" "ALX3_HUMAN.H11MO.0.D" ...
##   ..@ factors : list()
data("Probes.motif.hg38.450K")
dim(Probes.motif.hg38.450K)
## [1] 435391    771
str(Probes.motif.hg38.450K)
## Formal class 'ngCMatrix' [package "Matrix"] with 5 slots
##   ..@ i       : int [1:51916928] 0 5 6 9 12 13 14 17 18 20 ...
##   ..@ p       : int [1:772] 0 161625 175781 191186 206251 222385 265704 308216 348210 434158 ...
##   ..@ Dim     : int [1:2] 435391 771
##   ..@ Dimnames:List of 2
##   .. ..$ : chr [1:435391] "cg07568517" "cg17608500" "cg03125090" "cg26035915" ...
##   .. ..$ : chr [1:771] "AHR_HUMAN.H11MO.0.B" "AIRE_HUMAN.H11MO.0.C" "ALX1_HUMAN.H11MO.0.B" "ALX3_HUMAN.H11MO.0.D" ...
##   ..@ factors : list()
data("Probes.motif.hg19.EPIC")
dim(Probes.motif.hg19.EPIC)
## [1] 785511    771
str(Probes.motif.hg19.EPIC)
## Formal class 'ngCMatrix' [package "Matrix"] with 5 slots
##   ..@ i       : int [1:89283810] 0 1 5 14 15 17 18 24 28 32 ...
##   ..@ p       : int [1:772] 0 233791 266793 305520 342925 383080 483906 568097 653125 794659 ...
##   ..@ Dim     : int [1:2] 785511 771
##   ..@ Dimnames:List of 2
##   .. ..$ : chr [1:785511] "cg20457732" "cg21009965" "cg23154272" "cg00425213" ...
##   .. ..$ : chr [1:771] "AHR_HUMAN.H11MO.0.B" "AIRE_HUMAN.H11MO.0.C" "ALX1_HUMAN.H11MO.0.B" "ALX3_HUMAN.H11MO.0.D" ...
##   ..@ factors : list()
data("Probes.motif.hg38.EPIC")
dim(Probes.motif.hg38.EPIC)
## [1] 785511    771
str(Probes.motif.hg38.EPIC)
## Formal class 'ngCMatrix' [package "Matrix"] with 5 slots
##   ..@ i       : int [1:89270742] 0 1 3 6 8 11 15 16 17 20 ...
##   ..@ p       : int [1:772] 0 233791 266788 305500 342902 383049 483872 568051 653073 794331 ...
##   ..@ Dim     : int [1:2] 785511 771
##   ..@ Dimnames:List of 2
##   .. ..$ : chr [1:785511] "cg02230408" "cg00359465" "cg16743939" "cg20961824" ...
##   .. ..$ : chr [1:771] "AHR_HUMAN.H11MO.0.B" "AIRE_HUMAN.H11MO.0.C" "ALX1_HUMAN.H11MO.0.B" "ALX3_HUMAN.H11MO.0.D" ...
##   ..@ factors : list()

The following code was used to create the objects:

getInfiniumAnnotation <- function(plat = "450K", genome = "hg38"){
  message("Loading object: ",file)
  newenv <- new.env()
  if(plat == "EPIC" & genome == "hg19") data("EPIC.manifest", package = "ELMER.data",envir=newenv)
  if(plat == "EPIC" & genome == "hg38") data("EPIC.manifest.hg38", package = "ELMER.data",envir=newenv)
  if(plat == "450K" & genome == "hg19") data("hm450.manifest", package = "ELMER.data",envir=newenv)
  if(plat == "450K" & genome == "hg38") data("hm450.manifest.hg38", package = "ELMER.data",envir=newenv)
  annotation <- get(ls(newenv)[1],envir=newenv)   
  return(annotation)  
}
# To find for each probe the know motif we will use HOMER software (http://homer.salk.edu/homer/)
# Step:
# 1 - get DNA methylation probes annotation with the regions
# 2 - Make a bed file from it
# 3 - Execute section: Finding Instance of Specific Motifs from http://homer.salk.edu/homer/ngs/peakMotifs.html to the HOCOMOCO TF motifs
# Also, As HOMER is using more RAM than the available we will split the files in to 100k probes.
# Obs: for each probe we create a winddow of 500 bp (-size 500) around it. This might lead to false positives, but will not have false negatives.
# The false posives will be removed latter with some statistical tests.
TFBS.motif <- "http://hocomoco11.autosome.ru/final_bundle/hocomoco11/full/HUMAN/mono/HOCOMOCOv11_full_HUMAN_mono_homer_format_0.0001.motif"
if(!file.exists(basename(TFBS.motif))) downloader::download(TFBS.motif,basename(TFBS.motif))
for(plat in c("450K","EPIC")){
  for(gen in c("hg19","hg38")){
    
      file <- paste0(plat,gen,".txt")
      print(file)
      if(!file.exists(file)){
        # STEP 1
        gr <- getInfiniumAnnotation(plat = plat,genome =  gen)
        
        # This will remove masked probes. They have poor quality and might be arbitrarily positioned (Wanding Zhou)
        gr <- gr[!gr$MASK.general]

        df <- data.frame(seqnames=seqnames(gr),
                         starts=as.integer(start(gr)),
                         ends=end(gr),
                         names=names(gr),
                         scores=c(rep(".", length(gr))),
                         strands=strand(gr))
        step <- 10000 # nb of lines in each file. 10K was selected to not explode RAM
        n <- nrow(df)
        for(j in 0:floor(n/step)){
          # STEP 2
          file.aux <- paste0(plat,gen,"_",j,".bed")
          if(!file.exists(gsub(".bed",".txt",file.aux))){
            end <- ifelse(((j + 1) * step) > n, n,((j + 1) * step))
            write.table(df[((j * step) + 1):end,], file = file.aux, col.names = F, quote = F,row.names = F,sep = "\t")
            
            # STEP 3 use -mscore to get scores
            cmd <- paste("source ~/.bash_rc; annotatePeaks.pl" ,file.aux, gen, "-m", basename(TFBS.motif), "-size 500 -cpu 12 >", gsub(".bed",".txt",file.aux))
            system(cmd)
          }
        }
      }
      # We will merge the results from each file into one
      peaks <- NULL
      for(j in 0:floor(n/step)){
        aux <-  readr::read_tsv(paste0(plat,gen,"_",j,".txt"))
        colnames(aux)[1] <- "PeakID"
        if(is.null(peaks)) {
          peaks <- aux
        } else {
          peaks <- rbind(peaks, aux)
        }
      }
      readr::write_tsv(peaks,path=file,col_names = TRUE)
      print("DONE!")
      gc()
    }
  }
}

getMatrix <- function(filename) {
  motifs <- readr::read_tsv(file)
  # From 1 to 21 we have annotations then we have 640 motifs
  matrix <- Matrix::Matrix(0, nrow = nrow(motifs), ncol = ncol(motifs) - 21 ,sparse = TRUE)
  colnames(matrix) <- gsub(" Distance From Peak\\(sequence,strand,conservation\\)","",colnames(motifs)[-c(1:21)])
  rownames(matrix) <- motifs$PeakID
  matrix[!is.na(motifs[,-c(1:21)])] <- 1
  matrix <- as(matrix, "nsparseMatrix")
  return(matrix)
}

# This code will read the table with the motifs, save it as a sparce matrix
# and save all as a .rda that will be placed in ELMER.data
getMatrix <- function(filename) {
  motifs <- readr::read_tsv(file)
  # From 1 to 21 we have annotations then we have 640 motifs
  matrix <- Matrix::Matrix(0, nrow = nrow(motifs), ncol = 640,sparse = TRUE)
  colnames(matrix) <- gsub(" Distance From Peak\\(sequence,strand,conservation\\)","",colnames(motifs)[-c(1:21)])
  rownames(matrix) <- motifs$PeakID
  matrix[!is.na(motifs[,-c(1:21)])] <- 1
  matrix <- as(matrix, "nsparseMatrix")
  return(matrix)
}
for(plat in c("EPIC", "450K")){
  for(gen in c("hg19","hg38")){
    file <- paste0(plat,gen,".txt")
    if(file == "450Khg19.txt"){
      Probes.motif.hg19.450K <- getMatrix(file)
      save(Probes.motif.hg19.450K, file = "Probes.motif.hg19.450K.rda", compress = "xz")
      rm(Probes.motif.hg19.450K)
    } 
    if(file == "450Khg38.txt"){
      Probes.motif.hg38.450K <- getMatrix(file)
      save(Probes.motif.hg38.450K, file = "Probes.motif.hg38.450K.rda", compress = "xz")
      rm(Probes.motif.hg38.450K)
    }
    
    if(file == "EPIChg19.txt"){
      Probes.motif.hg19.EPIC <- getMatrix(file)
      save(Probes.motif.hg19.EPIC, file = "Probes.motif.hg19.EPIC.rda", compress = "xz")
      rm(Probes.motif.hg19.EPIC)
    }
    
    if(file == "EPIChg38.txt"){
      Probes.motif.hg38.EPIC <- getMatrix(file)
      save(Probes.motif.hg38.EPIC, file = "Probes.motif.hg38.EPIC.rda", compress = "xz")
      rm(Probes.motif.hg38.EPIC)
    }
  }
}
data("Probes.motif.hg19.450K")
as.data.frame(as.matrix(Probes.motif.hg19.450K[1:20,1:20])) %>% 
  datatable(options = list(scrollX = TRUE,pageLength = 5)) 
Homer versions
  • Software: v4.9.1
  • Genome hg19: v5.10
  • Genome hg38: v5.10

3 Session Information


sessionInfo()
## R version 3.4.3 (2017-11-30)
## 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=C                 
##  [9] LC_ADDRESS=C               LC_TELEPHONE=C            
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C       
## 
## attached base packages:
## [1] stats4    parallel  stats     graphics  grDevices utils     datasets 
## [8] methods   base     
## 
## other attached packages:
##  [1] Matrix_1.2-12       Biostrings_2.46.0   XVector_0.18.0     
##  [4] IRanges_2.12.0      S4Vectors_0.16.0    BiocGenerics_0.24.0
##  [7] dplyr_0.7.4         DT_0.2              ELMER.data_2.2.2   
## [10] BiocStyle_2.6.1    
## 
## loaded via a namespace (and not attached):
##  [1] Rcpp_0.12.14            compiler_3.4.3          GenomeInfoDb_1.14.0    
##  [4] bindr_0.1               bitops_1.0-6            tools_3.4.3            
##  [7] zlibbioc_1.24.0         digest_0.6.12           lattice_0.20-35        
## [10] jsonlite_1.5            evaluate_0.10.1         tibble_1.3.4           
## [13] pkgconfig_2.0.1         rlang_0.1.4             yaml_2.1.15            
## [16] bindrcpp_0.2            GenomeInfoDbData_0.99.1 stringr_1.2.0          
## [19] knitr_1.17              htmlwidgets_0.9         grid_3.4.3             
## [22] rprojroot_1.2           glue_1.2.0              R6_2.2.2               
## [25] rmarkdown_1.8           bookdown_0.5            magrittr_1.5           
## [28] backports_1.1.1         htmltools_0.3.6         GenomicRanges_1.30.0   
## [31] assertthat_0.2.0        stringi_1.1.6           RCurl_1.95-4.8

References

Heinz, Sven and Benner, Christopher and Spann, Nathanael and Bertolino, Eric and Lin, Yin C and Laslo, Peter and Cheng, Jason X and Murre, Cornelis and Singh, Harinder and Glass, Christopher K. 2010. “Simple Combinations of Lineage-Determining Transcription Factors Prime Cis-Regulatory Elements Required for Macrophage and B Cell Identities.”

Kulakovskiy, Ivan V and Medvedeva, Yulia A and Schaefer, Ulf and Kasianov, Artem S and Vorontsov, Ilya E and Bajic, Vladimir B and Makeev, Vsevolod. 2013. “HOCOMOCO a Comprehensive Collection of Human Transcription Factor Binding Sites Models.”

Wingender, E., Schoeps, T., Haubrock, M., & Dönitz, J. 2014. “TFClass a Classification of Human Transcription Factors and Their Rodent Orthologs.”

Yao, L., Shen, H., Laird, P. W., Farnham, P. J., & Berman, B. P. 2015. “Inferring Regulatory Element Landscapes and Transcription Factor Networks from Cancer Methylomes.”

Zhou, Wanding and Laird, Peter W and Shen, Hui. 2016. “Comprehensive Characterization, Annotation and Innovative Use of Infinium DNA Methylation BeadChip Probes.”